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