Arcane  4.1.12.0
Developer documentation
Loading...
Searching...
No Matches
VtkPolyhedralMeshIOService.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2026 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* VtkPolyhedralMeshIOService.cc (C) 2000-2025 */
9/* */
10/* Read/write fools for polyhedral mesh with vtk file format */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
14#include <iostream>
15#include <numeric>
16#include <functional>
17#include <memory>
18#include <iterator>
19
20#include <vtkUnstructuredGrid.h>
21#include <vtkUnstructuredGridReader.h>
22#include <vtkNew.h>
23#include <vtkCellIterator.h>
24#include <vtkIdTypeArray.h>
25#include <vtkCellData.h>
26#include <vtkPointData.h>
27#include <vtkDataSetAttributes.h>
28#include <vtkArrayDispatch.h>
29#include <vtkDataArrayAccessor.h>
30#include <vtkPolyDataReader.h>
31#include <vtkPolyData.h>
32#include <vtkSmartPointer.h>
33#include <vtkCellArray.h>
34#include <vtkVersion.h>
35#if VTK_VERSION_NUMBER >= 900000000
36using vtkIdType_generic = vtkIdType const;
37#else
38using vtkIdType_generic = vtkIdType;
39#endif
40
41#include <arccore/base/Ref.h>
42#include <arccore/base/String.h>
43#include <arccore/base/FatalErrorException.h>
44#include <vtkXMLUnstructuredGridReader.h>
45#include <vtkXMLPolyDataReader.h>
46
48#include "arcane/core/AbstractService.h"
49#include "arcane/core/ICaseMeshReader.h"
51#include "arcane/core/IMeshBuilder.h"
52#include "arcane/core/MeshBuildInfo.h"
53#include "arcane/core/IPrimaryMesh.h"
55#include "arcane/core/IMeshInitialAllocator.h"
56#include "arcane/core/IVariableMng.h"
58#include "arcane/utils/ITraceMng.h"
59#include "arcane/utils/UniqueArray.h"
60#include "arcane/utils/Real3.h"
61#include "arcane/mesh/CellFamily.h"
62#include "arcane/core/MeshVariableScalarRef.h"
63#include "arcane/core/MeshVariableArrayRef.h"
64
65#include "arcane/core/ItemAllocationInfo.h"
66#include "arcane/core/VariableBuildInfo.h"
67
68#include "arcane/utils/OStringStream.h"
69
70#include "arcane/core/IXmlDocumentHolder.h"
71#include "arcane/core/XmlNode.h"
72#include "arcane/core/internal/IVariableMngInternal.h"
73#include "arcane/core/datatype/DataTypeTraits.h"
74
75#include "arcane/std/VtkPolyhedralMeshIO_axl.h"
76
77/*---------------------------------------------------------------------------*/
78/*---------------------------------------------------------------------------*/
79
80namespace Arcane
81{
82
83/*---------------------------------------------------------------------------*/
84/*---------------------------------------------------------------------------*/
85
86namespace VtkPolyhedralTools
87{
89 {
90 bool failure = false;
91 String failure_message;
92 String info_message;
93 };
94
96 {
97 bool print_mesh_info = false;
98 bool print_debug_info = false;
99 };
100} // namespace VtkPolyhedralTools
101
102/*---------------------------------------------------------------------------*/
103/*---------------------------------------------------------------------------*/
104
105class VtkPolyhedralMeshIOService
106: public TraceAccessor
107{
108 public:
109
110 explicit VtkPolyhedralMeshIOService(ITraceMng* trace_mng, VtkPolyhedralTools::PrintInfoLevel print_info_level)
111 : TraceAccessor(trace_mng)
112 , m_print_info_level(print_info_level)
113 {}
114
115 class VtkReader
116 {
117
118 public:
119
120 explicit VtkReader(const String& filename, VtkPolyhedralTools::PrintInfoLevel print_info_level = VtkPolyhedralTools::PrintInfoLevel{});
121
122 VtkReader() = default;
123
124 static String supportedVtkExtensions() noexcept { return "vtk,vtu"; };
125
126 Int64ConstArrayView cellUids();
127 Int64ConstArrayView nodeUids();
128 Int64ConstArrayView faceUids();
129 Int64ConstArrayView edgeUids();
130
131 Integer nbNodes();
132
133 Int64ConstArrayView cellNodes();
134 Int32ConstArrayView cellNbNodes();
135
136 Int64ConstArrayView faceNodes();
137 Int32ConstArrayView faceNbNodes();
138
139 Int64ConstArrayView faceNodesInFaceMesh();
140 Int32ConstArrayView faceNbNodesInFaceMesh();
141
142 Int64ConstArrayView edgeNodes();
143 Int32ConstArrayView edgeNbNodes();
144
145 Int64ConstArrayView faceCells();
146 Int32ConstArrayView faceNbCells();
147
148 Int32ConstArrayView edgeNbCells();
149 Int64ConstArrayView edgeCells();
150
151 Int32ConstArrayView cellNbFaces();
152 Int64ConstArrayView cellFaces();
153
154 Int32ConstArrayView edgeNbFaces();
155 Int64ConstArrayView edgeFaces();
156
157 Int32ConstArrayView cellNbEdges();
158 Int64ConstArrayView cellEdges();
159
160 Int32ConstArrayView faceNbEdges();
161 Int64ConstArrayView faceEdges();
162
163 Int32ConstArrayView nodeNbCells();
164 Int64ConstArrayView nodeCells();
165
166 Int32ConstArrayView nodeNbFaces();
167 Int64ConstArrayView nodeFaces();
168
169 Int32ConstArrayView nodeNbEdges();
170 Int64ConstArrayView nodeEdges();
171
172 Real3ArrayView nodeCoords();
173
174 bool readHasFailed() const noexcept { return m_read_status.failure; }
175 VtkPolyhedralTools::ReadStatus readStatus() const noexcept { return m_read_status; }
176
177 vtkCellData* cellData();
178 vtkPointData* pointData();
179 vtkCellData* faceData();
180
181 bool isEmpty() const noexcept { return m_is_empty; }
182 bool doRead() const noexcept { return m_do_read; }
183
184 private:
185
186 bool m_is_empty = true;
187 bool m_do_read = false;
188 String m_filename;
189 VtkPolyhedralTools::PrintInfoLevel m_print_info_level;
190 VtkPolyhedralTools::ReadStatus m_read_status;
191 vtkNew<vtkUnstructuredGridReader> m_vtk_grid_reader;
192 vtkNew<vtkXMLUnstructuredGridReader> m_vtk_xml_grid_reader;
193 vtkNew<vtkPolyDataReader> m_vtk_face_grid_reader;
194 vtkNew<vtkXMLPolyDataReader> m_vtk_xml_face_grid_reader;
195 vtkUnstructuredGrid* m_vtk_grid = nullptr;
196 vtkPolyData* m_vtk_face_grid = nullptr;
197 Int64UniqueArray m_cell_uids, m_node_uids, m_face_uids, m_edge_uids;
198 Int64UniqueArray m_face_node_uids, m_edge_node_uids, m_cell_node_uids;
199 Int64UniqueArray m_face_cell_uids, m_edge_cell_uids, m_edge_face_uids;
200 Int64UniqueArray m_cell_face_uids, m_cell_edge_uids, m_face_edge_uids;
201 Int64UniqueArray m_node_cell_uids, m_node_face_uids, m_node_edge_uids;
202 Int32UniqueArray m_face_nb_nodes, m_cell_nb_nodes, m_face_nb_cells;
203 Int32UniqueArray m_edge_nb_cells, m_edge_nb_faces, m_cell_nb_faces;
204 Int32UniqueArray m_node_nb_cells, m_node_nb_faces, m_node_nb_edges;
205 Int32UniqueArray m_cell_nb_edges, m_face_nb_edges, m_face_uid_indexes;
206 Int32UniqueArray m_cell_face_indexes, m_edge_nb_nodes;
207 Int64UniqueArray m_face_node_uids_in_face_mesh;
208 Int32UniqueArray m_face_nb_nodes_in_face_mesh;
209 using NodeUidToIndexMap = Int32UniqueArray; // use a map when no longer possible to index with node uids
210 using FaceUidToIndexMap = Int32UniqueArray; // use a map when no longer possible to index with face uids
211 using EdgeUidToIndexMap = Int32UniqueArray; // use a map when no longer possible to index with edge uids
212 NodeUidToIndexMap m_node_uid_to_index;
213 Real3UniqueArray m_node_coordinates;
214 vtkCellData* m_cell_data = nullptr;
215 vtkPointData* m_point_data = nullptr;
216 vtkCellData* m_face_data = nullptr;
217 vtkCellArray* m_poly_data = nullptr; // to store faces from face mesh file
218
219 void _printMeshInfos() const;
220
221 std::pair<bool, Int32> _findFace(const Int64UniqueArray& sorted_face_nodes, const UniqueArray<Int64UniqueArray>& node_face_uids, const NodeUidToIndexMap& node_uid_to_index, const Int32UniqueArray& face_nb_nodes, const FaceUidToIndexMap& face_uid_to_index, const UniqueArray<Int32>& face_offsets, const Int64UniqueArray& face_node_uids);
222 template <typename Connectivity2DArray>
223 static void _flattenConnectivity(Connectivity2DArray connected_item_2darray, Int32Span nb_connected_item_per_source_item, Int64UniqueArray& connected_item_array);
224 void _readPlainTextVtkGrid(const String& filename);
225 void _readXlmVtkGrid(const String& filename);
226 void _checkVtkGrid() const;
227 void _readPlainTextVtkFaceGrid(const String& faces_filename);
228 void _readXmlVtkFaceGrid(const String& faces_filename);
229 void _readfaceNodesInFaceMesh();
230 };
231
232 public:
233
234 VtkPolyhedralTools::ReadStatus read(IPrimaryMesh* mesh, const String& filename, bool is_parallel_read)
235 {
237 // warning meaning of parallel_read not obvious :
238 // if is parallel_read => master reads + broadcast, if not all ranks read a pre-distributed mesh
239 bool do_read = is_parallel_read ? mesh->parallelMng()->isMasterIO() : true;
240 using VtkReaderPtr = std::unique_ptr<VtkReader>;
241 VtkReaderPtr reader = std::make_unique<VtkReader>();
242 if (do_read)
243 reader = std::make_unique<VtkReader>(filename, m_print_info_level);
244 if (reader->readHasFailed())
245 return reader->readStatus();
246 ItemAllocationInfo item_allocation_info;
247 _fillItemAllocationInfo(item_allocation_info, *reader);
248 auto polyhedral_mesh_allocator = mesh->initialAllocator()->polyhedralMeshAllocator();
249 polyhedral_mesh_allocator->allocateItems(item_allocation_info);
250 _readVariablesAndGroups(mesh, *reader);
251 return reader->readStatus();
252 }
253
254 private:
255
256 VtkPolyhedralTools::PrintInfoLevel m_print_info_level;
257
259 {
260 eDataType m_type = DT_Unknown;
261 bool is_array = false;
262 };
263
264 void _readVariablesAndGroups(IPrimaryMesh* mesh, VtkReader& reader);
265 void _createGroup(vtkDataArray* group_items, const String& group_name, IPrimaryMesh* mesh, IItemFamily* item_family, Int32ConstSpan vtkToArcaneLid) const;
266 VariableInfo _createVariable(vtkDataArray* item_values, const String& variable_name, IMesh* mesh, IItemFamily* item_family, Int32ConstSpan arcane_to_vtk_lids) const;
267
268 static void _fillItemAllocationInfo(ItemAllocationInfo& item_allocation_info, VtkReader& vtk_reader);
269 void _computeFaceVtkArcaneLidConversion(Int32Span face_vtk_to_arcane_lids, Int32Span arcane_to_vtk_lids, VtkPolyhedralMeshIOService::VtkReader& reader, IPrimaryMesh* mesh) const;
270 void _createEmptyVariablesAndGroups(IMesh* mesh, XmlNode::const_reference xml_node) const;
271 template <template <class> class VariableRootType, template <class> class ArrayVariableRootType>
272 void _createEmptyVariables(IMesh* mesh, const XmlNodeList& cell_variables_node, eItemKind item_kind) const;
273 void _createEmptyGroups(IMesh* mesh, const XmlNodeList& children, IItemFamily* item_family) const;
274};
275
276/*---------------------------------------------------------------------------*/
277/*---------------------------------------------------------------------------*/
278
279class VtkPolyhedralCaseMeshReader
281{
282 public:
283
284 class Builder : public IMeshBuilder
285 {
286 public:
287
288 explicit Builder(ITraceMng* tm, const CaseMeshReaderReadInfo& read_info, VtkPolyhedralTools::PrintInfoLevel print_info_level)
289 : m_trace_mng(tm)
290 , m_read_info(read_info)
291 , m_print_info_level(print_info_level)
292 {}
293
294 void fillMeshBuildInfo(MeshBuildInfo& build_info) override
295 {
296 build_info.addFactoryName("ArcanePolyhedralMeshFactory");
297 build_info.addNeedPartitioning(false);
298 MeshKind mk;
299 mk.setMeshStructure(eMeshStructure::Polyhedral);
300 build_info.addMeshKind(mk);
301 }
302
304 {
306 m_trace_mng->info() << "---Create Polyhedral mesh: " << pm->name() << "---";
307 m_trace_mng->info() << "--Read mesh file " << m_read_info.fileName();
308 VtkPolyhedralMeshIOService polyhedral_vtk_service{ m_trace_mng, m_print_info_level };
309 auto read_status = polyhedral_vtk_service.read(pm, m_read_info.fileName(), m_read_info.isParallelRead());
310 if (read_status.failure)
311 ARCANE_FATAL(read_status.failure_message);
312 m_trace_mng->info() << read_status.info_message;
313 }
314
315 private:
316
317 ITraceMng* m_trace_mng;
318 CaseMeshReaderReadInfo m_read_info;
319 VtkPolyhedralTools::PrintInfoLevel m_print_info_level;
320 };
321
324 {}
325
326 public:
327
329 {
330 IMeshBuilder* builder = nullptr;
331 if (VtkPolyhedralMeshIOService::VtkReader::supportedVtkExtensions().contains(read_info.format()))
332 builder = new Builder(traceMng(), read_info, VtkPolyhedralTools::PrintInfoLevel{ options()->getPrintMeshInfos(), options()->getPrintDebugInfos() });
333 return makeRef(builder);
334 }
335};
336
337/*---------------------------------------------------------------------------*/
338/*---------------------------------------------------------------------------*/
339
340ARCANE_REGISTER_SERVICE(VtkPolyhedralCaseMeshReader,
341 ServiceProperty("VtkPolyhedralCaseMeshReader", ST_CaseOption),
342 ARCANE_SERVICE_INTERFACE(ICaseMeshReader));
343
344/*---------------------------------------------------------------------------*/
345/*---------------------------------------------------------------------------*/
346
347void VtkPolyhedralMeshIOService::
348_fillItemAllocationInfo(ItemAllocationInfo& item_allocation_info, VtkReader& vtk_reader)
349{
350 auto nb_item_family = 4;
351 auto nb_connected_family = 3;
352 item_allocation_info.family_infos.resize(nb_item_family);
353 for (auto& family_info : item_allocation_info.family_infos) {
354 family_info.connected_family_infos.resize(nb_connected_family);
355 }
356 // Create regular item families and connectivities
357 auto& cell_family_info = item_allocation_info.family_infos[0];
358 cell_family_info.name = "Cell";
359 cell_family_info.item_kind = IK_Cell;
360 cell_family_info.item_uids = vtk_reader.cellUids();
361 auto& node_family_info = item_allocation_info.family_infos[1];
362 node_family_info.name = "Node";
363 node_family_info.item_kind = IK_Node;
364 node_family_info.item_uids = vtk_reader.nodeUids();
365 auto& face_family_info = item_allocation_info.family_infos[2];
366 face_family_info.name = "Face";
367 face_family_info.item_kind = IK_Face;
368 face_family_info.item_uids = vtk_reader.faceUids();
369 auto& edge_family_info = item_allocation_info.family_infos[3];
370 edge_family_info.name = "Edge";
371 edge_family_info.item_kind = IK_Edge;
372 edge_family_info.item_uids = vtk_reader.edgeUids();
373 // Cell to nodes connectivity
374 auto cell_connected_family_index = 0;
375 auto& cell_connected_node_family_info = cell_family_info.connected_family_infos[cell_connected_family_index++];
376 cell_connected_node_family_info.name = node_family_info.name;
377 cell_connected_node_family_info.item_kind = node_family_info.item_kind;
378 cell_connected_node_family_info.connectivity_name = "CellToNodes";
379 cell_connected_node_family_info.nb_connected_items_per_item = vtk_reader.cellNbNodes();
380 cell_connected_node_family_info.connected_items_uids = vtk_reader.cellNodes();
381 // Cell to faces connectivity
382 auto& cell_connected_face_family_info = cell_family_info.connected_family_infos[cell_connected_family_index++];
383 cell_connected_face_family_info.name = face_family_info.name;
384 cell_connected_face_family_info.item_kind = face_family_info.item_kind;
385 cell_connected_face_family_info.connectivity_name = "CellToFaces";
386 cell_connected_face_family_info.nb_connected_items_per_item = vtk_reader.cellNbFaces();
387 cell_connected_face_family_info.connected_items_uids = vtk_reader.cellFaces();
388 // Cell to edges connectivity
389 auto& cell_connected_edge_family_info = cell_family_info.connected_family_infos[cell_connected_family_index++];
390 cell_connected_edge_family_info.name = edge_family_info.name;
391 cell_connected_edge_family_info.item_kind = edge_family_info.item_kind;
392 cell_connected_edge_family_info.connectivity_name = "CellToEdges";
393 cell_connected_edge_family_info.nb_connected_items_per_item = vtk_reader.cellNbEdges();
394 cell_connected_edge_family_info.connected_items_uids = vtk_reader.cellEdges();
395 // Face to cells connectivity
396 auto face_connected_family_index = 0;
397 auto& face_connected_cell_family_info = face_family_info.connected_family_infos[face_connected_family_index++];
398 face_connected_cell_family_info.name = cell_family_info.name;
399 face_connected_cell_family_info.item_kind = cell_family_info.item_kind;
400 face_connected_cell_family_info.connectivity_name = "FaceToCells";
401 face_connected_cell_family_info.nb_connected_items_per_item = vtk_reader.faceNbCells();
402 face_connected_cell_family_info.connected_items_uids = vtk_reader.faceCells();
403 // Face to nodes connectivity
404 auto& face_connected_node_family_info = face_family_info.connected_family_infos[face_connected_family_index++];
405 face_connected_node_family_info.name = node_family_info.name;
406 face_connected_node_family_info.item_kind = node_family_info.item_kind;
407 face_connected_node_family_info.connectivity_name = "FaceToNodes";
408 face_connected_node_family_info.nb_connected_items_per_item = vtk_reader.faceNbNodes();
409 face_connected_node_family_info.connected_items_uids = vtk_reader.faceNodes();
410 // Face to edges connectivity
411 auto& face_connected_edge_family_info = face_family_info.connected_family_infos[face_connected_family_index];
412 face_connected_edge_family_info.name = edge_family_info.name;
413 face_connected_edge_family_info.item_kind = edge_family_info.item_kind;
414 face_connected_edge_family_info.connectivity_name = "FaceToEdges";
415 face_connected_edge_family_info.nb_connected_items_per_item = vtk_reader.faceNbEdges();
416 face_connected_edge_family_info.connected_items_uids = vtk_reader.faceEdges();
417 // Edge to cells connectivity
418 auto edge_connected_family_index = 0;
419 auto& edge_connected_cell_family_info = edge_family_info.connected_family_infos[edge_connected_family_index++];
420 edge_connected_cell_family_info.name = cell_family_info.name;
421 edge_connected_cell_family_info.item_kind = cell_family_info.item_kind;
422 edge_connected_cell_family_info.connectivity_name = "EdgeToCells";
423 edge_connected_cell_family_info.nb_connected_items_per_item = vtk_reader.edgeNbCells();
424 edge_connected_cell_family_info.connected_items_uids = vtk_reader.edgeCells();
425 // Edge to faces connectivity
426 auto& edge_connected_face_family_info = edge_family_info.connected_family_infos[edge_connected_family_index++];
427 edge_connected_face_family_info.name = face_family_info.name;
428 edge_connected_face_family_info.item_kind = face_family_info.item_kind;
429 edge_connected_face_family_info.connectivity_name = "EdgeToFaces";
430 edge_connected_face_family_info.nb_connected_items_per_item = vtk_reader.edgeNbFaces();
431 edge_connected_face_family_info.connected_items_uids = vtk_reader.edgeFaces();
432 // Edge to nodes connectivity
433 auto& edge_connected_node_family_info = edge_family_info.connected_family_infos[edge_connected_family_index++];
434 edge_connected_node_family_info.name = node_family_info.name;
435 edge_connected_node_family_info.item_kind = node_family_info.item_kind;
436 edge_connected_node_family_info.connectivity_name = "EdgeToNodes";
437 edge_connected_node_family_info.nb_connected_items_per_item = vtk_reader.edgeNbNodes();
438 edge_connected_node_family_info.connected_items_uids = vtk_reader.edgeNodes();
439 // Node to cells connectivity
440 auto node_connected_family_index = 0;
441 auto& node_connected_cell_family_info = node_family_info.connected_family_infos[node_connected_family_index++];
442 node_connected_cell_family_info.name = cell_family_info.name;
443 node_connected_cell_family_info.item_kind = cell_family_info.item_kind;
444 node_connected_cell_family_info.connectivity_name = "NodeToCells";
445 node_connected_cell_family_info.nb_connected_items_per_item = vtk_reader.nodeNbCells();
446 node_connected_cell_family_info.connected_items_uids = vtk_reader.nodeCells();
447 // Node to faces connectivity
448 auto& node_connected_face_family_info = node_family_info.connected_family_infos[node_connected_family_index++];
449 node_connected_face_family_info.name = face_family_info.name;
450 node_connected_face_family_info.item_kind = face_family_info.item_kind;
451 node_connected_face_family_info.connectivity_name = "NodeToFaces";
452 node_connected_face_family_info.nb_connected_items_per_item = vtk_reader.nodeNbFaces();
453 node_connected_face_family_info.connected_items_uids = vtk_reader.nodeFaces();
454 // Node to edges connectivity
455 auto& node_connected_edge_family_info = node_family_info.connected_family_infos[node_connected_family_index++];
456 node_connected_edge_family_info.name = edge_family_info.name;
457 node_connected_edge_family_info.item_kind = edge_family_info.item_kind;
458 node_connected_edge_family_info.connectivity_name = "NodeToEdges";
459 node_connected_edge_family_info.nb_connected_items_per_item = vtk_reader.nodeNbEdges();
460 node_connected_edge_family_info.connected_items_uids = vtk_reader.nodeEdges();
461 // Node coordinates
462 node_family_info.item_coordinates_variable_name = "NodeCoord";
463 node_family_info.item_coordinates = vtk_reader.nodeCoords();
464}
465
466/*---------------------------------------------------------------------------*/
467/*---------------------------------------------------------------------------*/
468
469void VtkPolyhedralMeshIOService::
470_readVariablesAndGroups(IPrimaryMesh* mesh, VtkReader& reader)
471{
472 // Register variable and group info in an xml for non reader ranks
473 OStringStream created_infos_str;
474 created_infos_str() << "<?xml version='1.0' ?>\n";
475 created_infos_str() << "<infos>";
476 // Cell data
477 if (auto* cell_data = reader.cellData(); cell_data) {
478 // Read cell groups and variables
479 Int32SharedArray vtk_to_arcane_lids(mesh->cellFamily()->nbItem());
480 std::iota(vtk_to_arcane_lids.begin(), vtk_to_arcane_lids.end(), 0);
481 Int32SharedArray arcane_to_vtk_lids{ vtk_to_arcane_lids };
482 for (auto array_index = 0; array_index < cell_data->GetNumberOfArrays(); ++array_index) {
483 auto* cell_array = cell_data->GetArray(array_index);
484 if (!cell_array)
485 continue;
486 if (String name = cell_array->GetName(); name.substring(0, 6) == "GROUP_") {
487 _createGroup(cell_array, name.substring(6), mesh, mesh->cellFamily(), vtk_to_arcane_lids.constSpan());
488 created_infos_str() << "<cell-group name='" << name.substring(6) << "'/>";
489 }
490 else {
491 auto var_info = _createVariable(cell_array, name, mesh, mesh->cellFamily(), arcane_to_vtk_lids);
492 created_infos_str() << "<cell-variable name='" << name << "' "
493 << " data-type='" << dataTypeName(var_info.m_type) << "' "
494 << " is_array='" << std::boolalpha << var_info.is_array << "'/>";
495 }
496 if (m_print_info_level.print_debug_info) {
497 debug(Trace::High) << "Reading property " << cell_array->GetName();
498 for (auto tuple_index = 0; tuple_index < cell_array->GetNumberOfTuples(); ++tuple_index) {
499 for (auto component_index = 0; component_index < cell_array->GetNumberOfComponents(); ++component_index) {
500 debug(Trace::High) << cell_array->GetName() << "[" << tuple_index << "][" << component_index << "] = " << cell_array->GetComponent(tuple_index, component_index);
501 }
502 }
503 }
504 }
505 }
506 // Node data
507 if (auto* point_data = reader.pointData(); point_data) {
508 // Read node groups and variables
509 Int32SharedArray vtk_to_arcane_lids(mesh->nodeFamily()->nbItem());
510 std::iota(vtk_to_arcane_lids.begin(), vtk_to_arcane_lids.end(), 0);
511 Int32SharedArray arcane_to_vtk_lids{ vtk_to_arcane_lids };
512 for (auto array_index = 0; array_index < point_data->GetNumberOfArrays(); ++array_index) {
513 auto* point_array = point_data->GetArray(array_index);
514 if (String name = point_array->GetName(); name.substring(0, 6) == "GROUP_") {
515 _createGroup(point_array, name.substring(6), mesh, mesh->nodeFamily(), vtk_to_arcane_lids.constSpan());
516 created_infos_str() << "<node-group name='" << name.substring(6) << "'/>";
517 }
518 else {
519 auto var_info = _createVariable(point_array, name, mesh, mesh->nodeFamily(), arcane_to_vtk_lids);
520 created_infos_str() << "<node-variable name='" << name << "' "
521 << " data-type='" << dataTypeName(var_info.m_type) << "' "
522 << " is_array='" << std::boolalpha << var_info.is_array << "'/>";
523 }
524 if (m_print_info_level.print_debug_info) {
525 debug(Trace::High) << "Reading property " << point_array->GetName();
526 for (auto tuple_index = 0; tuple_index < point_array->GetNumberOfTuples(); ++tuple_index) {
527 for (auto component_index = 0; component_index < point_array->GetNumberOfComponents(); ++component_index) {
528 debug(Trace::High) << point_array->GetName() << "[" << tuple_index << "][" << component_index << "] = " << point_array->GetComponent(tuple_index, component_index);
529 }
530 }
531 }
532 }
533 }
534 // Face data
535 if (auto* face_data = reader.faceData(); face_data) {
536 // Read face groups and variables
537 Int32UniqueArray vtk_to_Arcane_lids(mesh->faceFamily()->nbItem());
538 Int32UniqueArray arcane_to_vtk_lids(mesh->faceFamily()->nbItem());
539 _computeFaceVtkArcaneLidConversion(vtk_to_Arcane_lids, arcane_to_vtk_lids, reader, mesh);
540 for (auto array_index = 0; array_index < face_data->GetNumberOfArrays(); ++array_index) {
541 auto* face_array = face_data->GetArray(array_index);
542 if (String name = face_array->GetName(); name.substring(0, 6) == "GROUP_") {
543 _createGroup(face_array, name.substring(6), mesh, mesh->faceFamily(), vtk_to_Arcane_lids);
544 created_infos_str() << "<face-group name='" << name.substring(6) << "'/>";
545 }
546 else {
547 auto var_info = _createVariable(face_array, name, mesh, mesh->faceFamily(), arcane_to_vtk_lids);
548 created_infos_str() << "<face-variable name='" << name << "' "
549 << " data-type='" << dataTypeName(var_info.m_type) << "' "
550 << " is_array='" << std::boolalpha << var_info.is_array << "'/>";
551 }
552 if (m_print_info_level.print_debug_info) {
553 debug(Trace::High) << "Reading property " << face_array->GetName();
554 for (auto tuple_index = 0; tuple_index < face_array->GetNumberOfTuples(); ++tuple_index) {
555 for (auto component_index = 0; component_index < face_array->GetNumberOfComponents(); ++component_index) {
556 debug(Trace::High) << face_array->GetName() << "[" << tuple_index << "][" << component_index << "] = " << face_array->GetComponent(tuple_index, component_index);
557 }
558 }
559 }
560 }
561 }
562 created_infos_str() << "</infos>";
563 // Create empty group and variables in non reader subdomain
564 ByteUniqueArray bytes;
565 auto* pm = mesh->parallelMng();
566 if (!reader.isEmpty() && pm->isMasterIO()) {
567 String str = created_infos_str.str();
568 ByteConstArrayView bv = str.utf8();
569 Integer len = bv.size();
570 bytes.resize(len + 1);
571 bytes.copy(bv);
572 }
573 pm->broadcastMemoryBuffer(bytes, pm->masterIORank());
574 if (reader.isEmpty()) {
575 std::unique_ptr<IXmlDocumentHolder> doc(IXmlDocumentHolder::loadFromBuffer(bytes, "InternalBuffer", traceMng()));
576 XmlNode doc_node = doc->documentNode();
577 _createEmptyVariablesAndGroups(mesh, doc_node);
578 }
579}
580
581/*---------------------------------------------------------------------------*/
582/*---------------------------------------------------------------------------*/
583
584void VtkPolyhedralMeshIOService::
585_createGroup(vtkDataArray* group_items, const String& group_name, IPrimaryMesh* mesh, IItemFamily* item_family, Int32ConstSpan vtkToArcaneLid) const
586{
587 ARCANE_CHECK_POINTER(group_items);
589 ARCANE_CHECK_POINTER(item_family);
590 if (group_items->GetNumberOfComponents() != 1)
591 fatal() << String::format("Cannot create item group {0}. Group information in data property must be a scalar", group_name);
592 debug() << "Create group " << group_name;
593 Int32UniqueArray arcane_lids;
594 arcane_lids.reserve((int)group_items->GetNumberOfValues());
595 using GroupDispatcher = vtkArrayDispatch::DispatchByValueType<vtkArrayDispatch::Integrals>;
596 auto group_creator = [&arcane_lids, &vtkToArcaneLid](auto* array) {
597 vtkIdType numTuples = array->GetNumberOfTuples();
598 vtkDataArrayAccessor<std::remove_pointer_t<decltype(array)>> array_accessor{ array };
599 auto local_id = 0;
600 for (vtkIdType tupleIdx = 0; tupleIdx < numTuples; ++tupleIdx) {
601 auto value = array_accessor.Get(tupleIdx, 0);
602 if (value)
603 arcane_lids.push_back(vtkToArcaneLid[local_id]);
604 ++local_id;
605 }
606 };
607 if (!GroupDispatcher::Execute(group_items, group_creator))
608 ARCANE_FATAL("Cannot create item group {0}. Group information in data property must be an integral type", group_name);
609 debug() << " local ids for item group " << group_name << " " << arcane_lids; // to remove
610 item_family->createGroup(group_name, arcane_lids);
611}
612
613/*---------------------------------------------------------------------------*/
614/*---------------------------------------------------------------------------*/
615
616template <typename vtkType>
618{
619 using type = vtkType;
620};
621
622template <>
623struct ToArcaneType<float>
624{
625 using type = Real;
626};
627
628template <>
629struct ToArcaneType<long long>
630{
631 using type = Int64;
632};
633
634template <typename T> using to_arcane_type_t = typename ToArcaneType<T>::type;
635/*---------------------------------------------------------------------------*/
636
637VtkPolyhedralMeshIOService::VariableInfo VtkPolyhedralMeshIOService::
638_createVariable(vtkDataArray* item_values, const String& variable_name, IMesh* mesh, IItemFamily* item_family, Int32ConstSpan arcane_to_vtk_lids) const
639{
640 ARCANE_CHECK_POINTER(item_values);
642 ARCANE_CHECK_POINTER(item_family);
643 if (item_values->GetNumberOfTuples() != item_family->nbItem())
644 ARCANE_FATAL("Cannot create variable {0}, {1} values are given for {2} items in {3} family",
645 variable_name, item_values->GetNumberOfTuples(), item_family->nbItem(), item_family->name());
646 debug() << "Create mesh variable " << variable_name;
647 VariableInfo var_info;
648 auto variable_creator = [mesh, variable_name, item_family, arcane_to_vtk_lids, this, &var_info](auto* values) {
649 VariableBuildInfo vbi{ mesh, variable_name };
650 using ValueType = typename std::remove_pointer_t<decltype(values)>::ValueType;
651 auto* var = new ItemVariableScalarRefT<to_arcane_type_t<ValueType>>{ vbi, item_family->itemKind() };
653 vtkDataArrayAccessor<std::remove_pointer_t<decltype(values)>> values_accessor{ values };
654 ENUMERATE_ITEM (item, item_family->allItems()) {
655 (*var)[item] = (to_arcane_type_t<ValueType>)values_accessor.Get(arcane_to_vtk_lids[item.localId()], 0);
656 }
657 var_info.m_type = DataTypeTraitsT<to_arcane_type_t<ValueType>>::type();
658 var_info.is_array = false;
659 };
660 auto array_variable_creator = [mesh, variable_name, item_family, arcane_to_vtk_lids, this, &var_info](auto* values) {
661 VariableBuildInfo vbi{ mesh, variable_name };
662 using ValueType = typename std::remove_pointer_t<decltype(values)>::ValueType;
663 auto* var = new ItemVariableArrayRefT<to_arcane_type_t<ValueType>>{ vbi, item_family->itemKind() };
665 vtkDataArrayAccessor<std::remove_pointer_t<decltype(values)>> values_accessor{ values };
666 var->resize(values->GetNumberOfComponents());
667 ENUMERATE_ITEM (item, item_family->allItems()) {
668 auto index = 0;
669 for (auto& var_value : (*var)[item]) {
670 var_value = (to_arcane_type_t<ValueType>)values_accessor.Get(arcane_to_vtk_lids[item.localId()], index++);
671 }
672 }
673 var_info.m_type = DataTypeTraitsT<to_arcane_type_t<ValueType>>::type();
674 var_info.is_array = true;
675 };
676 // Restrict to int and real values
677 using ValueTypes = vtkTypeList_Create_6(double, float, int, long, long long, short);
678 using ArrayDispatcher = vtkArrayDispatch::DispatchByValueType<ValueTypes>;
679 // Create ScalarVariable
680 bool is_variable_created = false;
681 if (item_values->GetNumberOfComponents() == 1) {
682 is_variable_created = ArrayDispatcher::Execute(item_values, variable_creator);
683 }
684 // Create ArrayVariable
685 else { // ArrayVariable
686 is_variable_created = ArrayDispatcher::Execute(item_values, array_variable_creator);
687 }
688 if (!is_variable_created)
689 ARCANE_FATAL("Cannot create variable {0}, it's data type is not supported. Only real and integral types are supported", variable_name);
690 return var_info;
691}
692
693/*---------------------------------------------------------------------------*/
694/*---------------------------------------------------------------------------*/
695
696void VtkPolyhedralMeshIOService::
697_computeFaceVtkArcaneLidConversion(Int32Span face_vtk_to_arcane_lids, Int32Span arcane_to_vtk_lids, VtkPolyhedralMeshIOService::VtkReader& reader, IPrimaryMesh* mesh) const
698{
699 auto face_nodes_unique_ids = reader.faceNodesInFaceMesh();
700 auto face_nb_nodes = reader.faceNbNodesInFaceMesh();
701 auto current_face_index = 0;
702 auto current_face_index_in_face_nodes = 0;
703 Int32UniqueArray face_nodes_local_ids(face_nodes_unique_ids.size());
704 mesh->nodeFamily()->itemsUniqueIdToLocalId(face_nodes_local_ids, face_nodes_unique_ids);
705 for (auto current_face_nb_node : face_nb_nodes) {
706 auto current_face_nodes = face_nodes_local_ids.subConstView(current_face_index_in_face_nodes, current_face_nb_node);
707 current_face_index_in_face_nodes += current_face_nb_node;
708 Node face_first_node{ mesh->nodeFamily()->view()[current_face_nodes[0]] };
709 face_vtk_to_arcane_lids[current_face_index] = MeshUtils::getFaceFromNodesLocalId(face_first_node, current_face_nodes).localId();
710 ++current_face_index;
711 }
712 auto vtk_lid = 0;
713 for (auto arcane_lid : face_vtk_to_arcane_lids) {
714 arcane_to_vtk_lids[arcane_lid] = vtk_lid;
715 ++vtk_lid;
716 }
717}
718
719/*---------------------------------------------------------------------------*/
720/*---------------------------------------------------------------------------*/
721
722template <template <class> class VariableRootType>
723VariableRef* _createVar(IMesh* mesh, const String& var_name, const String& var_data_type_name, eItemKind var_kind)
724{
725 bool has_error = false;
726 eDataType var_data_type = dataTypeFromName(var_data_type_name.localstr(), has_error);
727 if (has_error) {
728 ARCANE_FATAL("Invalid data type name {0} for Variable creation in VtkPolyhedralMeshIOService");
729 }
730 switch (var_data_type) {
731 case DT_Int32:
732 return new VariableRootType<Int32>{ VariableBuildInfo{ mesh, var_name }, var_kind };
733 break;
734 case DT_Int64:
735 return new VariableRootType<Int64>{ VariableBuildInfo{ mesh, var_name }, var_kind };
736 break;
737 case DT_Real:
738 return new VariableRootType<Real>{ VariableBuildInfo{ mesh, var_name }, var_kind };
739 break;
740 default:
741 ARCANE_FATAL("Handle only DT_Int32, DT_Int64, DT_Real in VtkPolyhedralMeshIOService");
742 }
743}
744
745/*---------------------------------------------------------------------------*/
746/*---------------------------------------------------------------------------*/
747
748template <template <class> class VariableRootType, template <class> class ArrayVariableRootType>
749void VtkPolyhedralMeshIOService::
750_createEmptyVariables(IMesh* mesh, const XmlNodeList& item_variables_node, eItemKind item_kind) const
751{
752 ARCANE_CHECK_PTR(mesh);
753 {
754 for (XmlNode xnode : item_variables_node) {
755 String name = xnode.attrValue("name");
756 debug() << "Create mesh variable: " << name;
757 String data_type_name = xnode.attrValue("data-type");
758 bool is_array = xnode.attrValue("is_array") == "true";
759 VariableRef* var = nullptr;
760 if (is_array) {
761 var = _createVar<ArrayVariableRootType>(mesh, name, data_type_name, item_kind);
762 }
763 else {
764 var = _createVar<VariableRootType>(mesh, name, data_type_name, item_kind);
765 }
766 mesh->variableMng()->_internalApi()->addAutoDestroyVariable(var);
767 }
768 }
769}
770
771/*---------------------------------------------------------------------------*/
772/*---------------------------------------------------------------------------*/
773
774void VtkPolyhedralMeshIOService::
775_createEmptyGroups([[maybe_unused]] IMesh* mesh, const XmlNodeList& groups_node,
776 IItemFamily* item_family) const
777{
778 for (XmlNode xnode : groups_node) {
779 String name = xnode.attrValue("name");
780 info() << "Building group: " << name;
781 item_family->createGroup(name);
782 }
783}
784
785/*---------------------------------------------------------------------------*/
786/*---------------------------------------------------------------------------*/
787
788void VtkPolyhedralMeshIOService::
789_createEmptyVariablesAndGroups(IMesh* mesh, XmlNode::const_reference variable_and_group_info) const
790{
791 ARCANE_CHECK_PTR(mesh);
792 auto document_node = variable_and_group_info.documentElement();
793 _createEmptyVariables<ItemVariableScalarRefT, ItemVariableArrayRefT>(mesh, document_node.children("cell-variable"), IK_Cell);
794 _createEmptyVariables<ItemVariableScalarRefT, ItemVariableArrayRefT>(mesh, document_node.children("node-variable"), IK_Node);
795 _createEmptyVariables<ItemVariableScalarRefT, ItemVariableArrayRefT>(mesh, document_node.children("face-variable"), IK_Face);
796
797 _createEmptyGroups(mesh, document_node.children("cell-group"), mesh->itemFamily(IK_Cell));
798 _createEmptyGroups(mesh, document_node.children("node-group"), mesh->itemFamily(IK_Node));
799 _createEmptyGroups(mesh, document_node.children("face-group"), mesh->itemFamily(IK_Face));
800}
801
802/*---------------------------------------------------------------------------*/
803/*---------------------------------------------------------------------------*/
804
805VtkPolyhedralMeshIOService::VtkReader::
806VtkReader(const String& filename, VtkPolyhedralTools::PrintInfoLevel print_info_level)
807: m_filename{ filename }
808, m_print_info_level{ print_info_level }
809{
810 m_is_empty = false;
811 m_do_read = true;
812 if (filename.empty()) {
813 m_read_status.failure = true;
814 m_read_status.failure_message = "filename for polyhedral vtk mesh is empty.";
815 return;
816 }
817 if (filename.endsWith("vtk"))
818 _readPlainTextVtkGrid(filename);
819 else if (filename.endsWith("vtu"))
820 _readXlmVtkGrid(filename);
821 else {
822 m_read_status.failure = true;
823 m_read_status.failure_message = String::format("Unsupported vtk extension for file {0}. Supported vtk extension for Polyhedral meshes are {1}",
824 filename, VtkReader::supportedVtkExtensions());
825 }
826 // Check vtk grid exists and not empty
827 if (!m_vtk_grid) {
828 m_read_status.failure = true;
829 m_read_status.failure_message = String::format("Cannot read vtk polyhedral file {0}. Vtk grid was not created.", filename);
830 return;
831 }
832 if (m_vtk_grid->GetNumberOfCells() == 0) {
833 m_read_status.failure = true;
834 m_read_status.failure_message = String::format("Cannot read vtk polyhedral file {0}. No cells were found.", filename);
835 return;
836 }
837 if (!m_vtk_grid->GetFaces()) {
838 m_read_status.failure = true;
839 m_read_status.failure_message = String::format("The given mesh vtk file {0} is not a polyhedral mesh, cannot read it", filename);
840 return;
841 }
842
843 m_cell_data = m_vtk_grid->GetCellData();
844 m_point_data = m_vtk_grid->GetPointData();
845
846 // Read face info (for variables and groups) if present
847 String faces_filename = m_filename + "faces.vtk";
848 std::ifstream ifile(faces_filename.localstr());
849 if (ifile) {
850 _readPlainTextVtkFaceGrid(faces_filename);
851 }
852 else {
853 faces_filename = m_filename + "faces.vtp";
854 ifile = std::ifstream{ faces_filename.localstr() };
855 if (ifile) {
856 _readXmlVtkFaceGrid(faces_filename);
857 }
858 }
859
860 StringUniqueArray faces_filename_and_extension;
861 faces_filename.split(faces_filename_and_extension, '.');
862
863 if (!ifile) {
864 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.",
865 faces_filename_and_extension[0],
866 faces_filename_and_extension[1]);
867 }
868 else {
869 if (m_vtk_face_grid) { // Check face vtk grid exists and not empty
870 if (m_vtk_face_grid->GetNumberOfCells() == 0) {
871 m_read_status.failure = true;
872 m_read_status.failure_message = m_read_status.failure_message + String::format(" Error in reading face information for groups in mesh file {0} ", faces_filename);
873 }
874 else {
875 m_face_data = m_vtk_face_grid->GetCellData();
876 m_poly_data = m_vtk_face_grid->GetPolys();
877 }
878 }
879 else {
880 m_read_status.failure = true;
881 m_read_status.failure_message = m_read_status.failure_message + String::format("Face data could not be built from file {0}{1} (.vtk or .vtu).", faces_filename_and_extension[0], faces_filename_and_extension[1]);
882 }
883 }
884
885 if (m_print_info_level.print_mesh_info)
886 _printMeshInfos();
887}
888
889/*---------------------------------------------------------------------------*/
890/*---------------------------------------------------------------------------*/
891
892Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
893cellUids()
894{
895 if (!doRead())
896 return m_cell_uids;
897 if (m_cell_uids.empty()) {
898 _checkVtkGrid();
899 m_cell_uids.reserve(m_vtk_grid->GetNumberOfCells());
900 m_cell_nb_nodes.reserve(m_vtk_grid->GetNumberOfCells());
901 m_cell_node_uids.reserve(10 * m_vtk_grid->GetNumberOfCells()); // take a mean of 10 nodes per cell
902 auto* cell_iter = m_vtk_grid->NewCellIterator();
903 cell_iter->InitTraversal();
904 while (!cell_iter->IsDoneWithTraversal()) {
905 m_cell_uids.push_back(cell_iter->GetCellId());
906 m_cell_nb_nodes.push_back(Integer(cell_iter->GetNumberOfPoints()));
907 ArrayView<vtkIdType> cell_nodes{ Integer(cell_iter->GetNumberOfPoints()), cell_iter->GetPointIds()->GetPointer(0) };
908 std::for_each(cell_nodes.begin(), cell_nodes.end(), [this](auto uid) { this->m_cell_node_uids.push_back(uid); });
909 cell_iter->GoToNextCell();
910 }
911 }
912 return m_cell_uids;
913}
914
915/*---------------------------------------------------------------------------*/
916/*---------------------------------------------------------------------------*/
917
918Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
919nodeUids()
920{
921 if (!doRead())
922 return m_node_uids;
923 if (m_node_uids.empty()) {
924 _checkVtkGrid();
925 auto nb_nodes = m_vtk_grid->GetNumberOfPoints();
926 m_node_uids.resize(nb_nodes);
927 m_node_nb_cells.resize(nb_nodes);
928 m_node_cell_uids.reserve(8 * nb_nodes);
929 m_node_uid_to_index.resize(nb_nodes);
930 for (int node_index = 0; node_index < nb_nodes; ++node_index) {
931 Int64 node_uid = node_index;
932 m_node_uids[node_index] = node_uid;
933 auto cell_nodes = vtkIdList::New();
934 m_vtk_grid->GetPointCells(node_index, cell_nodes);
935 Int64Span cell_nodes_view((Int64*)cell_nodes->GetPointer(0), cell_nodes->GetNumberOfIds());
936 m_node_cell_uids.addRange(cell_nodes_view);
937 m_node_nb_cells[node_index] = (Int32)cell_nodes->GetNumberOfIds();
938 // uid and index might differ in the future (ex parallel read).
939 // This structure is created to be used in faceUids. Should be a map if node_uid is no longer the node index
940 m_node_uid_to_index[node_uid] = node_index;
941 }
942 }
943 return m_node_uids;
944}
945
946/*---------------------------------------------------------------------------*/
947/*---------------------------------------------------------------------------*/
948
949Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
950faceUids()
951{
952 // Needs nodeUids to be called first (no work done if already called)
953 nodeUids();
954
955 if (!doRead())
956 return m_face_uids;
957 if (!m_face_uids.empty())
958 return m_face_uids;
959 _checkVtkGrid();
960 auto* cell_iter = m_vtk_grid->NewCellIterator();
961 cell_iter->InitTraversal();
962 vtkIdType nb_face_estimation = 0;
963 while (!cell_iter->IsDoneWithTraversal()) {
964 vtkIdType cell_nb_faces = 0;
965 vtkIdType_generic* points{ nullptr };
966 m_vtk_grid->GetFaceStream(cell_iter->GetCellId(), cell_nb_faces, points);
967 nb_face_estimation += cell_nb_faces;
968 cell_iter->GoToNextCell();
969 }
970 m_face_uids.reserve(nb_face_estimation);
971 auto const* faces = m_vtk_grid->GetFaces();
972 // This array contains the face info per cells (cf. vtk file)
973 // first_cell_nb_faces first_cell_first_face_nb_nodes first_cell_first_face_node_1 ... first_cell_first_face_node_n first_cell_second_face_nb_nodes etc
974
975 if (!faces) {
976 ARCANE_FATAL("Mesh {0} is not polyhedral: faces are not defined", m_filename);
977 }
978 Int64 face_uid = 0;
979 auto face_info_size = faces->GetNumberOfValues();
980 m_face_node_uids.reserve(face_info_size);
981 m_face_nb_nodes.reserve(nb_face_estimation);
982 m_face_cell_uids.reserve(2 * nb_face_estimation);
983 m_face_nb_cells.reserve(nb_face_estimation);
984 m_cell_face_uids.reserve(8 * m_cell_uids.size()); // take a mean of 8 faces per cell
985 m_cell_nb_faces.resize(m_cell_uids.size(), 0);
986 m_cell_face_indexes.resize(m_cell_uids.size(), -1);
987 m_face_uid_indexes.resize(2 * nb_face_estimation, -1);
988 Int64UniqueArray current_face_nodes, sorted_current_face_nodes;
989 current_face_nodes.reserve(10);
990 sorted_current_face_nodes.reserve(10);
991 UniqueArray<UniqueArray<Int64>> node_faces(m_node_uids.size());
992 UniqueArray<Int32> face_offsets;
993 face_offsets.reserve(nb_face_estimation);
994 face_offsets.push_back(0);
995 FaceUidToIndexMap face_uid_to_index;
996 face_uid_to_index.reserve(nb_face_estimation);
997 auto cell_index = 0;
998 auto cell_face_index = 0;
999 auto global_face_index = 0;
1000 auto face_uid_index = 0;
1001 for (int face_info_index = 0; face_info_index < face_info_size; cell_index++) { // face data are given by cell
1002 auto current_cell_nb_faces = Int32(faces->GetValue(face_info_index++));
1003 m_cell_face_indexes[m_cell_uids[cell_index]] = cell_face_index;
1004 for (auto face_index = 0; face_index < current_cell_nb_faces; ++face_index, ++global_face_index) {
1005 auto current_face_nb_nodes = Int32(faces->GetValue(face_info_index++));
1006 m_cell_nb_faces[m_cell_uids[cell_index]] += 1;
1007 for (int node_index = 0; node_index < current_face_nb_nodes; ++node_index) {
1008 current_face_nodes.push_back(faces->GetValue(face_info_index++));
1009 }
1010 sorted_current_face_nodes.resize(current_face_nodes.size());
1011 auto is_front_cell = mesh_utils::reorderNodesOfFace(current_face_nodes, sorted_current_face_nodes);
1012 auto [is_face_found, existing_face_index] = _findFace(sorted_current_face_nodes, node_faces,
1013 m_node_uid_to_index, m_face_nb_nodes,
1014 face_uid_to_index, face_offsets, m_face_node_uids);
1015 if (!is_face_found) {
1016 for (auto node_uid : current_face_nodes) {
1017 node_faces[m_node_uid_to_index[node_uid]].push_back(face_uid);
1018 }
1019 m_cell_face_uids.push_back(face_uid);
1020 m_face_uids.push_back(face_uid); // todo parallel
1021 m_face_nb_nodes.push_back(current_face_nb_nodes);
1022 m_face_node_uids.addRange(sorted_current_face_nodes);
1023 m_face_nb_cells.push_back(1);
1024 m_face_uid_indexes[global_face_index] = face_uid_index;
1025 face_uid_to_index.push_back(face_uid_index);
1026 auto previous_offset = face_offsets.back();
1027 face_offsets.push_back(previous_offset + sorted_current_face_nodes.size());
1028 ++face_uid;
1029 ++face_uid_index;
1030 if (is_front_cell) {
1031 m_face_cell_uids.push_back(NULL_ITEM_UNIQUE_ID);
1032 m_face_cell_uids.push_back(m_cell_uids[cell_index]);
1033 }
1034 else {
1035 m_face_cell_uids.push_back(m_cell_uids[cell_index]);
1036 m_face_cell_uids.push_back(NULL_ITEM_UNIQUE_ID);
1037 }
1038 }
1039 else {
1040 m_cell_face_uids.push_back(m_face_uids[existing_face_index]);
1041 m_face_nb_cells[existing_face_index] += 1;
1042 m_face_uid_indexes[global_face_index] = existing_face_index;
1043 // add cell to face cell connectivity
1044 if (is_front_cell) {
1045 if (m_face_cell_uids[2 * existing_face_index + 1] != NULL_ITEM_UNIQUE_ID) {
1046 ARCANE_FATAL("Problem in face orientation, face uid {0}, nodes {1}, same orientation in cell {2} and {3}. Change mesh file.",
1047 m_face_uids[existing_face_index],
1048 current_face_nodes,
1049 m_face_cell_uids[2 * existing_face_index + 1],
1050 m_cell_uids[cell_index]);
1051 }
1052 m_face_cell_uids[2 * existing_face_index + 1] = m_cell_uids[cell_index];
1053 }
1054 else {
1055 if (m_face_cell_uids[2 * existing_face_index] != NULL_ITEM_UNIQUE_ID) {
1056 ARCANE_FATAL("Problem in face orientation, face uid {0}, nodes {1}, same orientation in cell {2} and {3}. Change mesh file.",
1057 m_face_uids[existing_face_index],
1058 current_face_nodes,
1059 m_face_cell_uids[2 * existing_face_index],
1060 m_cell_uids[cell_index]);
1061 }
1062 m_face_cell_uids[2 * existing_face_index] = m_cell_uids[cell_index];
1063 }
1064 }
1065 current_face_nodes.clear();
1066 sorted_current_face_nodes.clear();
1067 }
1068 cell_face_index += m_cell_nb_faces[m_cell_uids[cell_index]];
1069 }
1070 // fill node_face_uids and node_nb_faces from node_faces (array form [nb_nodes][nb_connected_faces])
1071 m_node_nb_faces.resize(m_node_uids.size(), 0);
1072 _flattenConnectivity(node_faces.constSpan(), m_node_nb_faces, m_node_face_uids);
1073
1074 if (m_print_info_level.print_debug_info) {
1075 std::cout << "================FACE NODES ==============" << std::endl;
1076 std::copy(m_face_node_uids.begin(), m_face_node_uids.end(), std::ostream_iterator<Int64>(std::cout, " "));
1077 std::cout << std::endl;
1078 std::copy(m_face_nb_nodes.begin(), m_face_nb_nodes.end(), std::ostream_iterator<Int64>(std::cout, " "));
1079 std::cout << std::endl;
1080 std::copy(m_cell_face_indexes.begin(), m_cell_face_indexes.end(), std::ostream_iterator<Int64>(std::cout, " "));
1081 std::cout << std::endl;
1082 }
1083
1084 return m_face_uids;
1085}
1086
1087/*---------------------------------------------------------------------------*/
1088/*---------------------------------------------------------------------------*/
1089
1090Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1091edgeUids()
1092{
1093 if (!doRead())
1094 return m_edge_uids;
1095 if (!m_edge_uids.empty())
1096 return m_edge_uids;
1097
1098 nodeUids(); // Needed to be called first. No work done if already called
1099
1100 _checkVtkGrid();
1101 m_edge_uids.reserve(2 * m_vtk_grid->GetNumberOfPoints());
1102 auto const* faces = m_vtk_grid->GetFaces();
1103 // This array contains the face info per cells (cf. vtk file)
1104 // first_cell_nb_faces first_cell_first_face_nb_nodes first_cell_first_face_node_1 ... first_cell_first_face_node_n first_cell_second_face_nb_nodes etc
1105
1106 if (!faces) {
1107 ARCANE_FATAL("Mesh {0} is not polyhedral: faces are not defined", m_filename);
1108 }
1109 Int64 edge_uid = 0;
1110 auto nb_edge_estimation = 2 * m_edge_uids.capacity();
1111 m_edge_node_uids.reserve(nb_edge_estimation);
1112 auto face_info_size = faces->GetNumberOfValues();
1113 auto cell_index = 0;
1114 auto global_face_index = 0;
1115 auto new_edge_index = 0;
1116 UniqueArray<std::set<Int64>> edge_cells;
1117 UniqueArray<Int64UniqueArray> edge_faces;
1118 edge_cells.reserve(m_edge_uids.capacity());
1119 edge_faces.reserve(m_edge_uids.capacity());
1120 m_cell_nb_edges.resize(m_cell_uids.size(), 0);
1121 m_cell_edge_uids.reserve(20 * m_cell_uids.size()); // choose a value of 20 edge per cell
1122 UniqueArray<std::set<Int64>> face_edges;
1123 face_edges.resize(m_face_uids.size());
1124 UniqueArray<std::set<Int64>> cell_edges;
1125 cell_edges.resize(m_cell_uids.size());
1126 UniqueArray<Int64UniqueArray> node_edges;
1127 node_edges.resize(m_node_uids.size());
1128 EdgeUidToIndexMap edge_uid_to_index;
1129 edge_uid_to_index.reserve(nb_edge_estimation);
1130 Int32UniqueArray edge_offsets;
1131 edge_offsets.reserve(nb_edge_estimation);
1132 edge_offsets.push_back(0);
1133 m_edge_nb_nodes.reserve(nb_edge_estimation);
1134 for (int face_info_index = 0; face_info_index < face_info_size; ++cell_index) {
1135 auto current_cell_nb_faces = Int32(faces->GetValue(face_info_index++));
1136 for (auto face_index = 0; face_index < current_cell_nb_faces; ++face_index, ++global_face_index) {
1137 auto current_face_nb_nodes = Int32(faces->GetValue(face_info_index++));
1138 auto first_face_node_uid = Int32(faces->GetValue(face_info_index));
1139 UniqueArray<Int64> current_edge(2), sorted_edge(2);
1140 for (int node_index = 0; node_index < current_face_nb_nodes - 1; ++node_index) {
1141 current_edge = UniqueArray<Int64>{ faces->GetValue(face_info_index++), faces->GetValue(face_info_index) };
1142 mesh_utils::reorderNodesOfFace(current_edge, sorted_edge); // works for edges
1143 auto [is_edge_found, existing_edge_index] = _findFace(sorted_edge, node_edges,
1144 m_node_uid_to_index,
1145 m_edge_nb_nodes,
1146 edge_uid_to_index, edge_offsets, m_edge_node_uids); // works for edges
1147 if (!is_edge_found) {
1148 m_cell_nb_edges[cell_index] += 1;
1149 face_edges[m_face_uid_indexes[global_face_index]].insert(edge_uid);
1150 cell_edges[cell_index].insert(edge_uid);
1151 for (auto node : current_edge) {
1152 node_edges[node].push_back(edge_uid);
1153 }
1154 edge_cells.push_back(std::set{ m_cell_uids[cell_index] });
1155 edge_faces.push_back(Int64UniqueArray{ m_cell_face_uids[m_cell_face_indexes[cell_index] + face_index] });
1156 m_edge_uids.push_back(edge_uid++); // todo parallel
1157 m_edge_node_uids.addRange(sorted_edge);
1158 edge_uid_to_index.push_back(new_edge_index);
1159 auto current_offset = edge_offsets.back();
1160 edge_offsets.push_back(current_offset + 2);
1161 m_edge_nb_nodes.push_back(2);
1162 ++new_edge_index;
1163 }
1164 else {
1165 edge_cells[existing_edge_index].insert(m_cell_uids[cell_index]);
1166 edge_faces[existing_edge_index].push_back(m_cell_face_uids[m_cell_face_indexes[cell_index] + face_index]);
1167 face_edges[m_face_uid_indexes[global_face_index]].insert(m_edge_uids[existing_edge_index]);
1168 cell_edges[cell_index].insert(m_edge_uids[existing_edge_index]);
1169 }
1170 }
1171 current_edge = UniqueArray<Int64>{ faces->GetValue(face_info_index++), first_face_node_uid };
1172 mesh_utils::reorderNodesOfFace(current_edge, sorted_edge); // works for edges
1173 auto [is_edge_found, existing_edge_index] = _findFace(sorted_edge, node_edges,
1174 m_node_uid_to_index,
1175 m_edge_nb_nodes,
1176 edge_uid_to_index, edge_offsets, m_edge_node_uids); // works for edges
1177 if (!is_edge_found) {
1178 m_cell_nb_edges[cell_index] += 1;
1179 edge_cells.push_back(std::set{ m_cell_uids[cell_index] });
1180 edge_faces.push_back(Int64UniqueArray{ m_cell_face_uids[m_cell_face_indexes[cell_index] + face_index] });
1181 face_edges[m_face_uid_indexes[global_face_index]].insert(edge_uid);
1182 cell_edges[cell_index].insert(edge_uid);
1183 for (auto node : current_edge) {
1184 node_edges[node].push_back(edge_uid);
1185 }
1186 m_edge_uids.push_back(edge_uid++); // todo parallel
1187 m_edge_node_uids.addRange(sorted_edge);
1188 edge_uid_to_index.push_back(new_edge_index);
1189 auto current_offset = edge_offsets.back();
1190 edge_offsets.push_back(current_offset + 2);
1191 m_edge_nb_nodes.push_back(2);
1192 ++new_edge_index;
1193 }
1194 else {
1195 edge_cells[existing_edge_index].insert(m_cell_uids[cell_index]);
1196 edge_faces[existing_edge_index].push_back(m_cell_face_uids[m_cell_face_indexes[cell_index] + face_index]);
1197 face_edges[m_face_uid_indexes[global_face_index]].insert(m_edge_uids[existing_edge_index]);
1198 cell_edges[cell_index].insert(m_edge_uids[existing_edge_index]);
1199 }
1200 }
1201 }
1202 // fill edge_cell_uids and edge_nb_cells from edge_cells (array form [nb_edges][nb_connected_cells])
1203 m_edge_nb_cells.resize(m_edge_uids.size(), 0);
1204 _flattenConnectivity(edge_cells.constSpan(), m_edge_nb_cells, m_edge_cell_uids);
1205
1206 // fill edge faces uids
1207 m_edge_nb_faces.resize(m_edge_uids.size(), 0);
1208 _flattenConnectivity(edge_faces.constSpan(), m_edge_nb_faces, m_edge_face_uids);
1209
1210 // fill face edge uids
1211 m_face_nb_edges.resize(m_face_uids.size(), 0);
1212 _flattenConnectivity(face_edges.constSpan(), m_face_nb_edges, m_face_edge_uids);
1213
1214 // fill cell edge uids
1215 m_cell_nb_edges.resize(m_cell_uids.size(), 0);
1216 _flattenConnectivity(cell_edges, m_cell_nb_edges, m_cell_edge_uids);
1217
1218 // fill node edge uids
1219 m_node_nb_edges.resize(m_node_uids.size(), 0);
1220 _flattenConnectivity(node_edges, m_node_nb_edges, m_node_edge_uids);
1221
1222 if (m_print_info_level.print_debug_info) {
1223 std::cout << "================EDGE NODES ==============" << std::endl;
1224 std::copy(m_edge_node_uids.begin(), m_edge_node_uids.end(), std::ostream_iterator<Int64>(std::cout, " "));
1225 std::cout << std::endl;
1226 std::cout << "================FACE EDGES ==============" << std::endl;
1227 std::copy(m_face_nb_edges.begin(), m_face_nb_edges.end(), std::ostream_iterator<Int32>(std::cout, " "));
1228 std::cout << std::endl;
1229 std::copy(m_face_edge_uids.begin(), m_face_edge_uids.end(), std::ostream_iterator<Int64>(std::cout, " "));
1230 std::cout << std::endl;
1231 std::cout << "================CELL EDGES ==============" << std::endl;
1232 std::copy(m_cell_nb_edges.begin(), m_cell_nb_edges.end(), std::ostream_iterator<Int32>(std::cout, " "));
1233 std::cout << std::endl;
1234 std::copy(m_cell_edge_uids.begin(), m_cell_edge_uids.end(), std::ostream_iterator<Int64>(std::cout, " "));
1235 std::cout << std::endl;
1236 }
1237 return m_edge_uids;
1238}
1239
1240/*---------------------------------------------------------------------------*/
1241/*---------------------------------------------------------------------------*/
1242
1243std::pair<bool, Int32> VtkPolyhedralMeshIOService::VtkReader::
1244_findFace(const Int64UniqueArray& sorted_face_nodes,
1245 const UniqueArray<Int64UniqueArray>& node_face_uids,
1246 const NodeUidToIndexMap& node_uid_to_index,
1247 const Int32UniqueArray& face_nb_nodes,
1248 const FaceUidToIndexMap& face_uid_to_index,
1249 const UniqueArray<Int32>& face_offsets,
1250 const Int64UniqueArray& face_node_uids)
1251{
1252 auto first_node_uid = sorted_face_nodes[0];
1253 auto first_node_index = node_uid_to_index[first_node_uid];
1254 // If the face already exists it has already been registered in node_face connectivity
1255 for (auto face_uid : node_face_uids[first_node_index]) {
1256 auto face_index = face_uid_to_index[face_uid];
1257 auto face_offset = face_offsets[face_index];
1258 auto face_nb_node = face_nb_nodes[face_index];
1259 if (face_nb_node == sorted_face_nodes.size()) {
1260 bool is_same_face = true;
1261 for (auto index = 0; index < face_nb_node; ++index) {
1262 if (sorted_face_nodes[index] != face_node_uids[face_offset + index]) {
1263 is_same_face = false;
1264 }
1265 }
1266 if (is_same_face)
1267 return { true, face_index };
1268 }
1269 }
1270 return { false, -1 };
1271}
1272
1273/*---------------------------------------------------------------------------*/
1274/*---------------------------------------------------------------------------*/
1275
1276Integer VtkPolyhedralMeshIOService::VtkReader::
1277nbNodes()
1278{
1279 if (m_node_uids.empty())
1280 nodeUids();
1281 return m_node_uids.size();
1282}
1283
1284/*---------------------------------------------------------------------------*/
1285/*---------------------------------------------------------------------------*/
1286
1287Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1288cellNodes()
1289{
1290 if (m_cell_node_uids.empty())
1291 cellUids();
1292 return m_cell_node_uids;
1293}
1294
1295/*---------------------------------------------------------------------------*/
1296/*---------------------------------------------------------------------------*/
1297
1298Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1299cellNbNodes()
1300{
1301 if (m_cell_nb_nodes.empty())
1302 cellUids();
1303 return m_cell_nb_nodes;
1304}
1305
1306/*---------------------------------------------------------------------------*/
1307/*---------------------------------------------------------------------------*/
1308
1309Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1310faceNodes()
1311{
1312 if (m_face_node_uids.empty())
1313 faceUids();
1314 return m_face_node_uids;
1315}
1316
1317/*---------------------------------------------------------------------------*/
1318/*---------------------------------------------------------------------------*/
1319
1320Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1321faceNbNodes()
1322{
1323 if (m_face_nb_nodes.empty())
1324 faceUids();
1325 return m_face_nb_nodes;
1326}
1327
1328/*---------------------------------------------------------------------------*/
1329/*---------------------------------------------------------------------------*/
1330
1331Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1332faceNodesInFaceMesh()
1333{
1334 if (m_face_node_uids_in_face_mesh.empty())
1335 _readfaceNodesInFaceMesh();
1336 return m_face_node_uids_in_face_mesh;
1337}
1338
1339/*---------------------------------------------------------------------------*/
1340/*---------------------------------------------------------------------------*/
1341
1342Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1343faceNbNodesInFaceMesh()
1344{
1345 if (m_face_nb_nodes_in_face_mesh.empty())
1346 _readfaceNodesInFaceMesh();
1347 return m_face_nb_nodes_in_face_mesh;
1348}
1349
1350/*---------------------------------------------------------------------------*/
1351/*---------------------------------------------------------------------------*/
1352
1353Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1354edgeNbNodes()
1355{
1356 if (m_edge_node_uids.empty())
1357 edgeUids();
1358 return m_edge_nb_nodes;
1359}
1360
1361/*---------------------------------------------------------------------------*/
1362/*---------------------------------------------------------------------------*/
1363
1364Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1365edgeNodes()
1366{
1367 if (m_edge_node_uids.empty())
1368 edgeUids();
1369 return m_edge_node_uids;
1370}
1371
1372/*---------------------------------------------------------------------------*/
1373/*---------------------------------------------------------------------------*/
1374
1375Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1376faceCells()
1377{
1378 if (m_face_cell_uids.empty())
1379 faceUids();
1380 // debug
1381 if (m_print_info_level.print_debug_info) {
1382 std::cout << "=================FACE CELLS================="
1383 << "\n";
1384 std::copy(m_face_cell_uids.begin(), m_face_cell_uids.end(), std::ostream_iterator<Int64>(std::cout, " "));
1385 std::cout << "\n";
1386 std::cout << "=================END FACE CELLS================="
1387 << "\n";
1388 }
1389 return m_face_cell_uids;
1390}
1391
1392/*---------------------------------------------------------------------------*/
1393/*---------------------------------------------------------------------------*/
1394
1395Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1396faceNbCells()
1397{
1398 if (m_face_nb_cells.empty())
1399 faceUids();
1400 return m_face_nb_cells;
1401}
1402
1403/*---------------------------------------------------------------------------*/
1404/*---------------------------------------------------------------------------*/
1405
1406Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1407edgeNbCells()
1408{
1409 if (m_edge_nb_cells.empty())
1410 edgeUids();
1411 return m_edge_nb_cells;
1412}
1413
1414/*---------------------------------------------------------------------------*/
1415/*---------------------------------------------------------------------------*/
1416
1417Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1418edgeCells()
1419{
1420 if (m_edge_cell_uids.empty())
1421 edgeUids();
1422 return m_edge_cell_uids;
1423}
1424
1425/*---------------------------------------------------------------------------*/
1426/*---------------------------------------------------------------------------*/
1427
1428Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1429cellNbFaces()
1430{
1431 if (m_cell_nb_faces.empty())
1432 faceUids();
1433 return m_cell_nb_faces;
1434}
1435
1436/*---------------------------------------------------------------------------*/
1437/*---------------------------------------------------------------------------*/
1438
1439Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1440cellFaces()
1441{
1442 if (m_cell_face_uids.empty())
1443 faceUids();
1444 return m_cell_face_uids;
1445}
1446
1447/*---------------------------------------------------------------------------*/
1448/*---------------------------------------------------------------------------*/
1449
1450Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1451edgeNbFaces()
1452{
1453 if (m_edge_nb_faces.empty())
1454 edgeUids();
1455 return m_edge_nb_faces;
1456}
1457
1458/*---------------------------------------------------------------------------*/
1459/*---------------------------------------------------------------------------*/
1460
1461Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1462edgeFaces()
1463{
1464 if (m_edge_face_uids.empty())
1465 edgeUids();
1466 return m_edge_face_uids;
1467}
1468
1469/*---------------------------------------------------------------------------*/
1470/*---------------------------------------------------------------------------*/
1471
1472Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1473cellNbEdges()
1474{
1475 if (m_cell_nb_edges.empty())
1476 edgeUids();
1477 return m_cell_nb_edges;
1478}
1479
1480/*---------------------------------------------------------------------------*/
1481/*---------------------------------------------------------------------------*/
1482
1483Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1484cellEdges()
1485{
1486 if (m_cell_edge_uids.empty())
1487 edgeUids();
1488 return m_cell_edge_uids;
1489}
1490
1491/*---------------------------------------------------------------------------*/
1492/*---------------------------------------------------------------------------*/
1493
1494Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1495faceNbEdges()
1496{
1497 if (m_face_nb_edges.empty())
1498 edgeUids();
1499 return m_face_nb_edges;
1500}
1501
1502/*---------------------------------------------------------------------------*/
1503/*---------------------------------------------------------------------------*/
1504
1505Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1506faceEdges()
1507{
1508 if (m_face_edge_uids.empty())
1509 edgeUids();
1510 return m_face_edge_uids;
1511}
1512
1513/*---------------------------------------------------------------------------*/
1514/*---------------------------------------------------------------------------*/
1515
1516template <typename Connectivity2DArray>
1517void VtkPolyhedralMeshIOService::VtkReader::
1518_flattenConnectivity(Connectivity2DArray connected_item_2darray,
1519 Int32Span nb_connected_item_per_source_item,
1520 Int64UniqueArray& connected_item_array)
1521{
1522 // fill nb_connected_item_per_source_items
1523 std::transform(connected_item_2darray.begin(), connected_item_2darray.end(), nb_connected_item_per_source_item.begin(), [](auto const& connected_items) {
1524 return connected_items.size();
1525 });
1526 // fill connected_item_array
1527 connected_item_array.reserve(std::accumulate(nb_connected_item_per_source_item.begin(), nb_connected_item_per_source_item.end(), 0));
1528 std::for_each(connected_item_2darray.begin(), connected_item_2darray.end(), [&connected_item_array](auto const& connected_items) {
1529 for (auto const& connected_item : connected_items) {
1530 connected_item_array.push_back(connected_item);
1531 }
1532 });
1533}
1534
1535/*---------------------------------------------------------------------------*/
1536/*---------------------------------------------------------------------------*/
1537
1538Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1539nodeNbCells()
1540{
1541 if (m_node_nb_cells.empty())
1542 nodeUids();
1543 return m_node_nb_cells;
1544}
1545
1546/*---------------------------------------------------------------------------*/
1547/*---------------------------------------------------------------------------*/
1548
1549Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1550nodeCells()
1551{
1552 if (m_node_cell_uids.empty())
1553 nodeUids();
1554 return m_node_cell_uids;
1555}
1556
1557/*---------------------------------------------------------------------------*/
1558/*---------------------------------------------------------------------------*/
1559
1560Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1561nodeNbFaces()
1562{
1563 if (m_node_nb_faces.empty())
1564 faceUids();
1565 return m_node_nb_faces;
1566}
1567
1568/*---------------------------------------------------------------------------*/
1569/*---------------------------------------------------------------------------*/
1570
1571Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1572nodeFaces()
1573{
1574 if (m_node_face_uids.empty())
1575 faceUids();
1576 return m_node_face_uids;
1577}
1578
1579/*---------------------------------------------------------------------------*/
1580/*---------------------------------------------------------------------------*/
1581
1582Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1583nodeNbEdges()
1584{
1585 if (m_node_nb_edges.empty())
1586 edgeUids();
1587 return m_node_nb_edges;
1588}
1589
1590/*---------------------------------------------------------------------------*/
1591/*---------------------------------------------------------------------------*/
1592
1593Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1594nodeEdges()
1595{
1596 if (m_node_edge_uids.empty())
1597 edgeUids();
1598 return m_node_edge_uids;
1599}
1600
1601/*---------------------------------------------------------------------------*/
1602/*---------------------------------------------------------------------------*/
1603
1604Real3ArrayView VtkPolyhedralMeshIOService::VtkReader::
1605nodeCoords()
1606{
1607 if (!doRead())
1608 return m_node_coordinates;
1609 if (m_node_coordinates.empty()) {
1610 _checkVtkGrid();
1611 auto point_coords = m_vtk_grid->GetPoints()->GetData();
1612 if (m_print_info_level.print_debug_info) {
1613 std::cout << "======= Point COORDS ====" << std::endl;
1614 std::ostringstream oss;
1615 point_coords->PrintSelf(oss, vtkIndent{ 2 });
1616 std::cout << oss.str() << std::endl;
1617 }
1618 auto nb_nodes = m_vtk_grid->GetNumberOfPoints();
1619 for (int i = 0; i < nb_nodes; ++i) {
1620 if (m_print_info_level.print_debug_info) {
1621 std::cout << "==========current point coordinates : ( ";
1622 std::cout << *(point_coords->GetTuple(i)) << " , ";
1623 std::cout << *(point_coords->GetTuple(i) + 1) << " , ";
1624 std::cout << *(point_coords->GetTuple(i) + 2) << " ) ===" << std::endl;
1625 }
1626 m_node_coordinates.add({ *(point_coords->GetTuple(i)),
1627 *(point_coords->GetTuple(i) + 1),
1628 *(point_coords->GetTuple(i) + 2) });
1629 }
1630 }
1631 return m_node_coordinates;
1632}
1633
1634/*---------------------------------------------------------------------------*/
1635/*---------------------------------------------------------------------------*/
1636
1637vtkCellData* VtkPolyhedralMeshIOService::VtkReader::
1638cellData()
1639{
1640 return m_cell_data;
1641}
1642
1643/*---------------------------------------------------------------------------*/
1644/*---------------------------------------------------------------------------*/
1645
1646vtkPointData* VtkPolyhedralMeshIOService::VtkReader::
1647pointData()
1648{
1649 return m_point_data;
1650}
1651
1652/*---------------------------------------------------------------------------*/
1653/*---------------------------------------------------------------------------*/
1654
1655void VtkPolyhedralMeshIOService::VtkReader::
1656_printMeshInfos() const
1657{
1658 _checkVtkGrid();
1659 std::cout << "-- VTK GRID READ "
1660 << " NB CELLS " << m_vtk_grid->GetNumberOfCells() << std::endl;
1661 // Parse cells
1662 auto* cell_iter = m_vtk_grid->vtkDataSet::NewCellIterator();
1663 cell_iter->InitTraversal();
1664 vtkIdType_generic* cell_faces{ nullptr };
1665 vtkIdType nb_faces = 0;
1666 while (!cell_iter->IsDoneWithTraversal()) {
1667 std::cout << "---- visiting cell id " << cell_iter->GetCellId() << std::endl;
1668 std::cout << "---- cell number of faces " << cell_iter->GetNumberOfFaces() << std::endl;
1669 std::cout << "---- cell number of points " << cell_iter->GetNumberOfPoints() << std::endl;
1670 m_vtk_grid->GetFaceStream(cell_iter->GetCellId(), nb_faces, cell_faces);
1671 for (auto iface = 0; iface < nb_faces; ++iface) {
1672 auto face_nb_nodes = *cell_faces++;
1673 std::cout << "---- has face with " << face_nb_nodes << " nodes. Node ids : ";
1674 for (int inode = 0; inode < face_nb_nodes; ++inode) {
1675 std::cout << *cell_faces++ << " ";
1676 }
1677 std::cout << std::endl;
1678 }
1679 cell_iter->GoToNextCell();
1680 }
1681}
1682
1683/*---------------------------------------------------------------------------*/
1684/*---------------------------------------------------------------------------*/
1685
1686vtkCellData* VtkPolyhedralMeshIOService::VtkReader::faceData()
1687{
1688 return m_face_data;
1689}
1690/*---------------------------------------------------------------------------*/
1691/*---------------------------------------------------------------------------*/
1692
1693/*---------------------------------------------------------------------------*/
1694/*---------------------------------------------------------------------------*/
1695
1696void VtkPolyhedralMeshIOService::VtkReader::
1697_readPlainTextVtkGrid(const String& filename)
1698{
1699 m_vtk_grid_reader->SetFileName(filename.localstr());
1700 m_vtk_grid_reader->ReadAllScalarsOn();
1701 m_vtk_grid_reader->Update();
1702 m_vtk_grid = m_vtk_grid_reader->GetOutput();
1703}
1704
1705/*---------------------------------------------------------------------------*/
1706/*---------------------------------------------------------------------------*/
1707
1708void VtkPolyhedralMeshIOService::VtkReader::
1709_readXlmVtkGrid(const String& filename)
1710{
1711 m_vtk_xml_grid_reader->SetFileName(filename.localstr());
1712 m_vtk_xml_grid_reader->Update();
1713 m_vtk_grid = m_vtk_xml_grid_reader->GetOutput();
1714}
1715
1716/*---------------------------------------------------------------------------*/
1717/*---------------------------------------------------------------------------*/
1718
1719void VtkPolyhedralMeshIOService::VtkReader::
1720_readPlainTextVtkFaceGrid(const String& faces_filename)
1721{
1722 m_vtk_face_grid_reader->SetFileName(faces_filename.localstr());
1723 m_vtk_face_grid_reader->ReadAllScalarsOn();
1724 m_vtk_face_grid_reader->Update();
1725 m_vtk_face_grid = m_vtk_face_grid_reader->GetOutput();
1726}
1727
1728/*---------------------------------------------------------------------------*/
1729/*---------------------------------------------------------------------------*/
1730
1731void VtkPolyhedralMeshIOService::VtkReader::
1732_readXmlVtkFaceGrid(const String& faces_filename)
1733{
1734 m_vtk_xml_face_grid_reader->SetFileName(faces_filename.localstr());
1735 m_vtk_xml_face_grid_reader->Update();
1736 m_vtk_face_grid = m_vtk_xml_face_grid_reader->GetOutput();
1737}
1738
1739/*---------------------------------------------------------------------------*/
1740/*---------------------------------------------------------------------------*/
1741
1742void VtkPolyhedralMeshIOService::VtkReader::
1743_readfaceNodesInFaceMesh()
1744{
1745 m_face_nb_nodes_in_face_mesh.resize(m_poly_data->GetNumberOfCells());
1746 m_face_node_uids_in_face_mesh.reserve(m_poly_data->GetNumberOfCells() * m_poly_data->GetMaxCellSize());
1747 m_poly_data->InitTraversal();
1748 vtkIdType face_nb_nodes;
1749 vtkIdType_generic* face_nodes;
1750
1751 auto face_nb_node_index = 0;
1752 Int64UniqueArray current_face_node_uids;
1753 Int64UniqueArray reordered_current_face_node_uids;
1754 current_face_node_uids.reserve(m_poly_data->GetMaxCellSize());
1755 reordered_current_face_node_uids.reserve(m_poly_data->GetMaxCellSize());
1756 Int64UniqueArray reordered_face_node_uids(m_poly_data->GetMaxCellSize());
1757 while (m_poly_data->GetNextCell(face_nb_nodes, face_nodes)) {
1758 m_face_nb_nodes_in_face_mesh[face_nb_node_index] = face_nb_nodes;
1759 ConstArrayView<vtkIdType> face_nodes_view(face_nb_nodes, face_nodes);
1760 current_face_node_uids.resize(face_nb_nodes);
1761 reordered_current_face_node_uids.resize(face_nb_nodes);
1762 std::copy(face_nodes_view.begin(), face_nodes_view.end(), current_face_node_uids.begin());
1763 MeshUtils::reorderNodesOfFace(current_face_node_uids, reordered_current_face_node_uids);
1764 std::copy(reordered_current_face_node_uids.begin(), reordered_current_face_node_uids.end(), std::back_inserter(m_face_node_uids_in_face_mesh));
1765 ++face_nb_node_index;
1766 }
1767}
1768/*---------------------------------------------------------------------------*/
1769/*---------------------------------------------------------------------------*/
1770
1771void VtkPolyhedralMeshIOService::VtkReader::
1772_checkVtkGrid() const
1773{
1774 if (!m_vtk_grid)
1775 ARCANE_FATAL("Polyhedral vtk grid not loaded. Cannot continue.");
1776}
1777
1778/*---------------------------------------------------------------------------*/
1779/*---------------------------------------------------------------------------*/
1780
1781} // End namespace Arcane
1782
1783/*---------------------------------------------------------------------------*/
1784/*---------------------------------------------------------------------------*/
Arcane configuration file.
#define ARCANE_CHECK_POINTER(ptr)
Macro returning the pointer ptr if it is not null or throwing an exception if it is null.
#define ARCANE_FATAL(...)
Macro throwing a FatalErrorException.
Declarations of Arcane's general types.
#define ENUMERATE_ITEM(name, group)
Generic enumerator for a node group.
Utility functions for the mesh.
This file contains the various service factories and macros for registering services.
#define ARCANE_SERVICE_INTERFACE(ainterface)
Macro to declare an interface when registering a service.
Management of references to a C++ class.
Generation de la classe de base du Service.
CaseOptionsVtkPolyhedralMeshIO * options() const
Options du jeu de données du service.
ArcaneVtkPolyhedralMeshIOObject(const Arcane::ServiceBuildInfo &sbi)
Constructeur.
void reserve(Int64 new_capacity)
Reserves memory for new_capacity elements.
Necessary information for reading a mesh file.
Interface of an entity family.
Definition IItemFamily.h:83
virtual String name() const =0
Family name.
virtual Integer nbItem() const =0
Number of entities.
virtual String name() const =0
Mesh name.
Interface of a mesh creation/reading service.
virtual IVariableMng * variableMng() const =0
Associated variable manager.
virtual void addAutoDestroyVariable(VariableRef *var)=0
Adds the variable to the list of variables that are kept until the end of execution.
virtual IVariableMngInternal * _internalApi()=0
Internal Arcane API.
static IXmlDocumentHolder * loadFromBuffer(Span< const Byte > buffer, const String &name, ITraceMng *tm)
Loads an XML document.
Definition DomUtils.cc:425
Parameters necessary for building a mesh.
MeshBuildInfo & addNeedPartitioning(bool v)
Indicates whether the generator needs to call a partitioner.
MeshBuildInfo & addMeshKind(const MeshKind &v)
Sets the mesh characteristics.
MeshBuildInfo & addFactoryName(const String &factory_name)
Sets the factory name to create this mesh.
Characteristics of a mesh.
Definition MeshKind.h:113
Reference to an instance.
Structure containing the information to create a service.
TraceAccessor(ITraceMng *m)
Constructs an accessor via the trace manager m.
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flow for a debug message.
TraceMessage fatal() const
Flow for a fatal error message.
TraceMessage info() const
Flow for an information message.
ITraceMng * traceMng() const
Trace manager.
1D data vector with value semantics (STL style).
Parameters necessary for building a variable.
Information characterizing a variable.
Reference to a variable.
Definition VariableRef.h:56
void fillMeshBuildInfo(MeshBuildInfo &build_info) override
Fills build_info with the necessary information to create the mesh.
void allocateMeshItems(IPrimaryMesh *pm) override
Allocates the mesh entities managed by this service.
Ref< IMeshBuilder > createBuilder(const CaseMeshReaderReadInfo &read_info) const override
Returns a builder to create and read the mesh whose information is specified in read_info.
List of nodes of a DOM tree.
Definition XmlNodeList.h:36
const value_type & const_reference
Type constant reference of an element in the array.
Definition XmlNode.h:67
#define ARCANE_REGISTER_SERVICE(aclass, a_service_property,...)
Macro for registering a service.
Integer len(const char *s)
Returns the length of the string s.
-- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature --
Span< Int64 > Int64Span
C equivalent of a 1D array of 64-bit integers.
Definition UtilsTypes.h:520
UniqueArray< Int64 > Int64UniqueArray
Dynamic 1D array of 64-bit integers.
Definition UtilsTypes.h:339
std::int64_t Int64
Signed integer type of 64 bits.
ArrayView< Real3 > Real3ArrayView
C equivalent of a 1D array of Real3.
Definition UtilsTypes.h:467
Int32 Integer
Type representing an integer.
UniqueArray< Real3 > Real3UniqueArray
Dynamic 1D array of rank 3 vectors.
Definition UtilsTypes.h:363
ConstArrayView< Int32 > Int32ConstArrayView
C equivalent of a 1D array of 32-bit integers.
Definition UtilsTypes.h:482
@ Polyhedral
Polyhedral mesh.
Definition MeshKind.h:38
@ ST_CaseOption
The service is used at the dataset level.
ConstArrayView< Int64 > Int64ConstArrayView
C equivalent of a 1D array of 64-bit integers.
Definition UtilsTypes.h:480
UniqueArray< Byte > ByteUniqueArray
Dynamic 1D array of characters.
Definition UtilsTypes.h:335
SharedArray< Int32 > Int32SharedArray
Dynamic 1D array of 32-bit integers.
Definition UtilsTypes.h:381
UniqueArray< Int32 > Int32UniqueArray
Dynamic 1D array of 32-bit integers.
Definition UtilsTypes.h:341
eItemKind
Mesh entity type.
@ IK_Node
Node mesh entity.
@ IK_Cell
Cell mesh entity.
@ IK_Face
Face mesh entity.
@ IK_Edge
Edge mesh entity.
double Real
Type representing a real number.
ConstArrayView< Byte > ByteConstArrayView
C equivalent of a 1D array of characters.
Definition UtilsTypes.h:476
auto makeRef(InstanceType *t) -> Ref< InstanceType >
Creates a reference on a pointer.
UniqueArray< String > StringUniqueArray
Dynamic 1D array of strings.
Definition UtilsTypes.h:359
Span< Int32 > Int32Span
C equivalent of a 1D array of 32-bit integers.
Definition UtilsTypes.h:522
eDataType
Data type.
Definition DataTypes.h:41
@ DT_Int32
32-bit integer data type
Definition DataTypes.h:45
@ DT_Int64
64-bit integer data type
Definition DataTypes.h:46
@ DT_Unknown
Unknown or uninitialized data type.
Definition DataTypes.h:58
@ DT_Real
Real data type.
Definition DataTypes.h:43
Span< const Int32 > Int32ConstSpan
Read-only view of a 1D array of 32-bit integers.
Definition UtilsTypes.h:554
ARCANE_DATATYPE_EXPORT eDataType dataTypeFromName(const char *name, bool &has_error)
Finds the type associated with name.
Definition DataTypes.cc:94
const char * dataTypeName(eDataType type)
Data type name.
Definition DataTypes.cc:72
std::int32_t Int32
Signed integer type of 32 bits.