Arcane  v3.14.10.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-2024 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-2024 */
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
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
33#include <arccore/base/Ref.h>
34#include <arccore/base/String.h>
35#include <arccore/base/FatalErrorException.h>
36#include <vtkXMLUnstructuredGridReader.h>
37#include <vtkXMLPolyDataReader.h>
38
40#include "arcane/core/AbstractService.h"
41#include "arcane/core/ICaseMeshReader.h"
43#include "arcane/core/IMeshBuilder.h"
44#include "arcane/core/MeshBuildInfo.h"
45#include "arcane/core/IPrimaryMesh.h"
47#include "arcane/core/IMeshInitialAllocator.h"
49#include "arcane/utils/ITraceMng.h"
50#include "arcane/utils/UniqueArray.h"
51#include "arcane/utils/Real3.h"
52#include "arcane/mesh/CellFamily.h"
53#include "arcane/core/MeshVariableScalarRef.h"
54#include "arcane/core/MeshVariableArrayRef.h"
55
56#include "arcane/core/ItemAllocationInfo.h"
57#include "arcane/core/VariableBuildInfo.h"
58
59#include "arcane/std/VtkPolyhedralMeshIO_axl.h"
60
61/*---------------------------------------------------------------------------*/
62/*---------------------------------------------------------------------------*/
63
64namespace Arcane
65{
66
67namespace VtkPolyhedralTools
68{
70 {
71 bool failure = false;
72 String failure_message;
73 String info_message;
74 };
75
77 {
78 bool print_mesh_info = false;
79 bool print_debug_info = false;
80 };
81} // namespace VtkPolyhedralTools
82
84: public TraceAccessor
85{
86 public:
87
89 : TraceAccessor(trace_mng)
90 , m_print_info_level(print_info_level)
91 {}
92
94 {
95
96 public:
97
99
100 static String supportedVtkExtensions() noexcept {return "vtk,vtu";};
101
102 Int64ConstArrayView cellUids();
103 Int64ConstArrayView nodeUids();
104 Int64ConstArrayView faceUids();
105 Int64ConstArrayView edgeUids();
106
107 Integer nbNodes();
108
109 Int64ConstArrayView cellNodes();
110 Int32ConstArrayView cellNbNodes();
111
112 Int64ConstArrayView faceNodes();
113 Int32ConstArrayView faceNbNodes();
114
115 Int64ConstArrayView edgeNodes();
116 Int32ConstArrayView edgeNbNodes();
117
118 Int64ConstArrayView faceCells();
119 Int32ConstArrayView faceNbCells();
120
121 Int32ConstArrayView edgeNbCells();
122 Int64ConstArrayView edgeCells();
123
124 Int32ConstArrayView cellNbFaces();
125 Int64ConstArrayView cellFaces();
126
127 Int32ConstArrayView edgeNbFaces();
128 Int64ConstArrayView edgeFaces();
129
130 Int32ConstArrayView cellNbEdges();
131 Int64ConstArrayView cellEdges();
132
133 Int32ConstArrayView faceNbEdges();
134 Int64ConstArrayView faceEdges();
135
136 Int32ConstArrayView nodeNbCells();
137 Int64ConstArrayView nodeCells();
138
139 Int32ConstArrayView nodeNbFaces();
140 Int64ConstArrayView nodeFaces();
141
142 Int32ConstArrayView nodeNbEdges();
143 Int64ConstArrayView nodeEdges();
144
145 Real3ArrayView nodeCoords();
146
147 bool readHasFailed() const noexcept { return m_read_status.failure; }
148 VtkPolyhedralTools::ReadStatus readStatus() const noexcept { return m_read_status; }
149
150 vtkCellData* cellData();
151 vtkPointData* pointData();
152 vtkCellData* faceData();
153
154 private:
155
156 const String& m_filename;
157 VtkPolyhedralTools::PrintInfoLevel m_print_info_level;
158 VtkPolyhedralTools::ReadStatus m_read_status;
159 vtkNew<vtkUnstructuredGridReader> m_vtk_grid_reader;
160 vtkNew<vtkXMLUnstructuredGridReader> m_vtk_xml_grid_reader;
161 vtkNew<vtkPolyDataReader> m_vtk_face_grid_reader;
162 vtkNew<vtkXMLPolyDataReader> m_vtk_xml_face_grid_reader;
163 vtkUnstructuredGrid* m_vtk_grid = nullptr;
164 vtkPolyData* m_vtk_face_grid = nullptr;
165 Int64UniqueArray m_cell_uids, m_node_uids, m_face_uids, m_edge_uids;
166 Int64UniqueArray m_face_node_uids, m_edge_node_uids, m_cell_node_uids;
167 Int64UniqueArray m_face_cell_uids, m_edge_cell_uids, m_edge_face_uids;
168 Int64UniqueArray m_cell_face_uids, m_cell_edge_uids, m_face_edge_uids;
169 Int64UniqueArray m_node_cell_uids, m_node_face_uids, m_node_edge_uids;
170 Int32UniqueArray m_face_nb_nodes, m_cell_nb_nodes, m_face_nb_cells;
171 Int32UniqueArray m_edge_nb_cells, m_edge_nb_faces, m_cell_nb_faces;
172 Int32UniqueArray m_node_nb_cells, m_node_nb_faces, m_node_nb_edges;
173 Int32UniqueArray m_cell_nb_edges, m_face_nb_edges, m_face_uid_indexes;
174 Int32UniqueArray m_cell_face_indexes, m_edge_nb_nodes;
175 using NodeUidToIndexMap = Int32UniqueArray; // use a map when no longer possible to index with node uids
176 using FaceUidToIndexMap = Int32UniqueArray; // use a map when no longer possible to index with face uids
177 using EdgeUidToIndexMap = Int32UniqueArray; // use a map when no longer possible to index with edge uids
178 NodeUidToIndexMap m_node_uid_to_index;
179 Real3UniqueArray m_node_coordinates;
180 vtkCellData* m_cell_data = nullptr;
181 vtkPointData* m_point_data = nullptr;
182 vtkCellData* m_face_data = nullptr;
183
184 void _printMeshInfos() const;
185
187 template <typename Connectivity2DArray>
189 void _readPlainTextVtkGrid(const String& filename);
190 void _readXlmVtkGrid(const String& filename);
191 void _checkVtkGrid() const;
192 void _readPlainTextVtkFaceGrid(const String& faces_filename);
193 void _readXmlVtkFaceGrid(const String& faces_filename);
194 };
195
196 public:
197
198 VtkPolyhedralTools::ReadStatus read(IPrimaryMesh* mesh, const String& filename)
199 {
201 VtkReader reader{ filename, m_print_info_level };
202 if (reader.readHasFailed())
203 return reader.readStatus();
205 _fillItemAllocationInfo(item_allocation_info, reader);
206 auto polyhedral_mesh_allocator = mesh->initialAllocator()->polyhedralMeshAllocator();
208 _readVariablesAndGroups(mesh, reader);
209 return reader.readStatus();
210 }
211
212 private:
213
214 UniqueArray<VariableRef*> m_read_variables;
215 VtkPolyhedralTools::PrintInfoLevel m_print_info_level;
216
217 void _readVariablesAndGroups(IPrimaryMesh* mesh, VtkReader& reader);
218 void _createGroup(vtkDataArray* group_items, const String& group_name, IPrimaryMesh* mesh, IItemFamily* item_family, Int32ConstSpan vtkToArcaneLid) const;
219 void _createVariable(vtkDataArray* item_values, const String& variable_name, IMesh* mesh, IItemFamily* item_family, Int32ConstSpan arcane_to_vtk_lids);
220
221 static void _fillItemAllocationInfo(ItemAllocationInfo& item_allocation_info, VtkReader& vtk_reader);
222 void _computeFaceVtkArcaneLidConversion(Int32Span face_vtk_to_arcane_lids, Int32Span arcane_to_vtk_lids, VtkPolyhedralMeshIOService::VtkReader& reader, IPrimaryMesh* mesh) const;
223};
224
225/*---------------------------------------------------------------------------*/
226/*---------------------------------------------------------------------------*/
227
230{
231 public:
232
233 class Builder : public IMeshBuilder
234 {
235 public:
236
238 : m_trace_mng(tm)
239 , m_read_info(read_info)
240 , m_print_info_level(print_info_level)
241 {}
242
244 {
245 build_info.addFactoryName("ArcanePolyhedralMeshFactory");
246 build_info.addNeedPartitioning(false);
247 MeshKind mk;
248 mk.setMeshStructure(eMeshStructure::Polyhedral);
249 build_info.addMeshKind(mk);
250 }
251
253 {
255 m_trace_mng->info() << "---Create Polyhedral mesh: " << pm->name() << "---";
256 m_trace_mng->info() << "--Read mesh file " << m_read_info.fileName();
257 VtkPolyhedralMeshIOService polyhedral_vtk_service{ m_trace_mng, m_print_info_level };
258 auto read_status = polyhedral_vtk_service.read(pm, m_read_info.fileName());
259 if (read_status.failure)
260 ARCANE_FATAL(read_status.failure_message);
261 m_trace_mng->info() << read_status.info_message;
262 }
263
264 private:
265
266 ITraceMng* m_trace_mng;
267 CaseMeshReaderReadInfo m_read_info;
268 VtkPolyhedralTools::PrintInfoLevel m_print_info_level;
269 };
270
273 {}
274
275 public:
276
278 {
279 IMeshBuilder* builder = nullptr;
280 if (VtkPolyhedralMeshIOService::VtkReader::supportedVtkExtensions().contains(read_info.format()))
281 builder = new Builder(traceMng(), read_info, VtkPolyhedralTools::PrintInfoLevel{ options()->getPrintMeshInfos(), options()->getPrintDebugInfos() });
282 return makeRef(builder);
283 }
284};
285
286/*---------------------------------------------------------------------------*/
287/*---------------------------------------------------------------------------*/
288
289ARCANE_REGISTER_SERVICE(VtkPolyhedralCaseMeshReader,
290 ServiceProperty("VtkPolyhedralCaseMeshReader", ST_CaseOption),
291 ARCANE_SERVICE_INTERFACE(ICaseMeshReader));
292
293/*---------------------------------------------------------------------------*/
294/*---------------------------------------------------------------------------*/
295
296void VtkPolyhedralMeshIOService::
297_fillItemAllocationInfo(ItemAllocationInfo& item_allocation_info, VtkReader& vtk_reader)
298{
299 auto nb_item_family = 4;
300 auto nb_connected_family = 3;
301 item_allocation_info.family_infos.resize(nb_item_family);
302 for (auto& family_info : item_allocation_info.family_infos) {
303 family_info.connected_family_info.resize(nb_connected_family);
304 }
305 // Create regular item families and connectivities
306 auto& cell_family_info = item_allocation_info.family_infos[0];
307 cell_family_info.name = "Cell";
308 cell_family_info.item_kind = IK_Cell;
309 cell_family_info.item_uids = vtk_reader.cellUids();
310 auto& node_family_info = item_allocation_info.family_infos[1];
311 node_family_info.name = "Node";
312 node_family_info.item_kind = IK_Node;
313 node_family_info.item_uids = vtk_reader.nodeUids();
314 auto& face_family_info = item_allocation_info.family_infos[2];
315 face_family_info.name = "Face";
316 face_family_info.item_kind = IK_Face;
317 face_family_info.item_uids = vtk_reader.faceUids();
318 auto& edge_family_info = item_allocation_info.family_infos[3];
319 edge_family_info.name = "Edge";
320 edge_family_info.item_kind = IK_Edge;
321 edge_family_info.item_uids = vtk_reader.edgeUids();
322 // Cell to nodes connectivity
323 auto cell_connected_family_index = 0;
324 auto& cell_connected_node_family_info = cell_family_info.connected_family_info[cell_connected_family_index++];
325 cell_connected_node_family_info.name = node_family_info.name;
326 cell_connected_node_family_info.item_kind = node_family_info.item_kind;
327 cell_connected_node_family_info.connectivity_name = "CellToNodes";
328 cell_connected_node_family_info.nb_connected_items_per_item = vtk_reader.cellNbNodes();
329 cell_connected_node_family_info.connected_items_uids = vtk_reader.cellNodes();
330 // Cell to faces connectivity
331 auto& cell_connected_face_family_info = cell_family_info.connected_family_info[cell_connected_family_index++];
332 cell_connected_face_family_info.name = face_family_info.name;
333 cell_connected_face_family_info.item_kind = face_family_info.item_kind;
334 cell_connected_face_family_info.connectivity_name = "CellToFaces";
335 cell_connected_face_family_info.nb_connected_items_per_item = vtk_reader.cellNbFaces();
336 cell_connected_face_family_info.connected_items_uids = vtk_reader.cellFaces();
337 // Cell to edges connectivity
338 auto& cell_connected_edge_family_info = cell_family_info.connected_family_info[cell_connected_family_index++];
339 cell_connected_edge_family_info.name = edge_family_info.name;
340 cell_connected_edge_family_info.item_kind = edge_family_info.item_kind;
341 cell_connected_edge_family_info.connectivity_name = "CellToEdges";
342 cell_connected_edge_family_info.nb_connected_items_per_item = vtk_reader.cellNbEdges();
343 cell_connected_edge_family_info.connected_items_uids = vtk_reader.cellEdges();
344 // Face to cells connectivity
345 auto face_connected_family_index = 0;
346 auto& face_connected_cell_family_info = face_family_info.connected_family_info[face_connected_family_index++];
347 face_connected_cell_family_info.name = cell_family_info.name;
348 face_connected_cell_family_info.item_kind = cell_family_info.item_kind;
349 face_connected_cell_family_info.connectivity_name = "FaceToCells";
350 face_connected_cell_family_info.nb_connected_items_per_item = vtk_reader.faceNbCells();
351 face_connected_cell_family_info.connected_items_uids = vtk_reader.faceCells();
352 // Face to nodes connectivity
353 auto& face_connected_node_family_info = face_family_info.connected_family_info[face_connected_family_index++];
354 face_connected_node_family_info.name = node_family_info.name;
355 face_connected_node_family_info.item_kind = node_family_info.item_kind;
356 face_connected_node_family_info.connectivity_name = "FaceToNodes";
357 face_connected_node_family_info.nb_connected_items_per_item = vtk_reader.faceNbNodes();
358 face_connected_node_family_info.connected_items_uids = vtk_reader.faceNodes();
359 // Face to edges connectivity
360 auto& face_connected_edge_family_info = face_family_info.connected_family_info[face_connected_family_index];
361 face_connected_edge_family_info.name = edge_family_info.name;
362 face_connected_edge_family_info.item_kind = edge_family_info.item_kind;
363 face_connected_edge_family_info.connectivity_name = "FaceToEdges";
364 face_connected_edge_family_info.nb_connected_items_per_item = vtk_reader.faceNbEdges();
365 face_connected_edge_family_info.connected_items_uids = vtk_reader.faceEdges();
366 // Edge to cells connectivity
367 auto edge_connected_family_index = 0;
368 auto& edge_connected_cell_family_info = edge_family_info.connected_family_info[edge_connected_family_index++];
369 edge_connected_cell_family_info.name = cell_family_info.name;
370 edge_connected_cell_family_info.item_kind = cell_family_info.item_kind;
371 edge_connected_cell_family_info.connectivity_name = "EdgeToCells";
372 edge_connected_cell_family_info.nb_connected_items_per_item = vtk_reader.edgeNbCells();
373 edge_connected_cell_family_info.connected_items_uids = vtk_reader.edgeCells();
374 // Edge to faces connectivity
375 auto& edge_connected_face_family_info = edge_family_info.connected_family_info[edge_connected_family_index++];
376 edge_connected_face_family_info.name = face_family_info.name;
377 edge_connected_face_family_info.item_kind = face_family_info.item_kind;
378 edge_connected_face_family_info.connectivity_name = "EdgeToFaces";
379 edge_connected_face_family_info.nb_connected_items_per_item = vtk_reader.edgeNbFaces();
380 edge_connected_face_family_info.connected_items_uids = vtk_reader.edgeFaces();
381 // Edge to nodes connectivity
382 auto& edge_connected_node_family_info = edge_family_info.connected_family_info[edge_connected_family_index++];
383 edge_connected_node_family_info.name = node_family_info.name;
384 edge_connected_node_family_info.item_kind = node_family_info.item_kind;
385 edge_connected_node_family_info.connectivity_name = "EdgeToNodes";
386 edge_connected_node_family_info.nb_connected_items_per_item = vtk_reader.edgeNbNodes();
387 edge_connected_node_family_info.connected_items_uids = vtk_reader.edgeNodes();
388 // Node to cells connectivity
389 auto node_connected_family_index = 0;
390 auto& node_connected_cell_family_info = node_family_info.connected_family_info[node_connected_family_index++];
391 node_connected_cell_family_info.name = cell_family_info.name;
392 node_connected_cell_family_info.item_kind = cell_family_info.item_kind;
393 node_connected_cell_family_info.connectivity_name = "NodeToCells";
394 node_connected_cell_family_info.nb_connected_items_per_item = vtk_reader.nodeNbCells();
395 node_connected_cell_family_info.connected_items_uids = vtk_reader.nodeCells();
396 // Node to faces connectivity
397 auto& node_connected_face_family_info = node_family_info.connected_family_info[node_connected_family_index++];
398 node_connected_face_family_info.name = face_family_info.name;
399 node_connected_face_family_info.item_kind = face_family_info.item_kind;
400 node_connected_face_family_info.connectivity_name = "NodeToFaces";
401 node_connected_face_family_info.nb_connected_items_per_item = vtk_reader.nodeNbFaces();
402 node_connected_face_family_info.connected_items_uids = vtk_reader.nodeFaces();
403 // Node to edges connectivity
404 auto& node_connected_edge_family_info = node_family_info.connected_family_info[node_connected_family_index++];
405 node_connected_edge_family_info.name = edge_family_info.name;
406 node_connected_edge_family_info.item_kind = edge_family_info.item_kind;
407 node_connected_edge_family_info.connectivity_name = "NodeToEdges";
408 node_connected_edge_family_info.nb_connected_items_per_item = vtk_reader.nodeNbEdges();
409 node_connected_edge_family_info.connected_items_uids = vtk_reader.nodeEdges();
410 // Node coordinates
411 node_family_info.item_coordinates_variable_name = "NodeCoord";
412 node_family_info.item_coordinates = vtk_reader.nodeCoords();
413}
414
415/*---------------------------------------------------------------------------*/
416/*---------------------------------------------------------------------------*/
417
418void VtkPolyhedralMeshIOService::
419_readVariablesAndGroups(IPrimaryMesh* mesh, VtkReader& reader)
420{
421 // Cell data
422 if (auto* cell_data = reader.cellData(); cell_data) {
423 // Read cell groups and variables
424 Int32SharedArray vtk_to_arcane_lids(mesh->cellFamily()->nbItem());
425 std::iota(vtk_to_arcane_lids.begin(), vtk_to_arcane_lids.end(), 0);
426 Int32SharedArray arcane_to_vtk_lids{ vtk_to_arcane_lids };
427 for (auto array_index = 0; array_index < cell_data->GetNumberOfArrays(); ++array_index) {
428 auto* cell_array = cell_data->GetArray(array_index);
429 if (!cell_array)
430 continue;
431 if (String name = cell_array->GetName(); name.substring(0, 6) == "GROUP_")
432 _createGroup(cell_array, name.substring(6), mesh, mesh->cellFamily(), vtk_to_arcane_lids.constSpan());
433 else
434 _createVariable(cell_array, name, mesh, mesh->cellFamily(), arcane_to_vtk_lids);
435 if (m_print_info_level.print_debug_info) {
436 debug(Trace::High) << "Reading property " << cell_array->GetName();
437 for (auto tuple_index = 0; tuple_index < cell_array->GetNumberOfTuples(); ++tuple_index) {
438 for (auto component_index = 0; component_index < cell_array->GetNumberOfComponents(); ++component_index) {
439 debug(Trace::High) << cell_array->GetName() << "[" << tuple_index << "][" << component_index << "] = " << cell_array->GetComponent(tuple_index, component_index);
440 }
441 }
442 }
443 }
444 }
445 // Node data
446 if (auto* point_data = reader.pointData(); point_data) {
447 // Read node groups and variables
448 Int32SharedArray vtk_to_arcane_lids(mesh->nodeFamily()->nbItem());
449 std::iota(vtk_to_arcane_lids.begin(), vtk_to_arcane_lids.end(), 0);
450 Int32SharedArray arcane_to_vtk_lids{ vtk_to_arcane_lids };
451 for (auto array_index = 0; array_index < point_data->GetNumberOfArrays(); ++array_index) {
452 auto* point_array = point_data->GetArray(array_index);
453 if (String name = point_array->GetName(); name.substring(0, 6) == "GROUP_")
454 _createGroup(point_array, name.substring(6), mesh, mesh->nodeFamily(), vtk_to_arcane_lids.constSpan());
455 else
456 _createVariable(point_array, name, mesh, mesh->nodeFamily(), arcane_to_vtk_lids);
457 if (m_print_info_level.print_debug_info) {
458 debug(Trace::High) << "Reading property " << point_array->GetName();
459 for (auto tuple_index = 0; tuple_index < point_array->GetNumberOfTuples(); ++tuple_index) {
460 for (auto component_index = 0; component_index < point_array->GetNumberOfComponents(); ++component_index) {
461 debug(Trace::High) << point_array->GetName() << "[" << tuple_index << "][" << component_index << "] = " << point_array->GetComponent(tuple_index, component_index);
462 }
463 }
464 }
465 }
466 }
467 // Face data
468 if (auto* face_data = reader.faceData(); face_data) {
469 // Read face groups and variables
470 Int32UniqueArray vtk_to_Arcane_lids(mesh->faceFamily()->nbItem());
471 Int32UniqueArray arcane_to_vtk_lids(mesh->faceFamily()->nbItem());
472 _computeFaceVtkArcaneLidConversion(vtk_to_Arcane_lids, arcane_to_vtk_lids, reader, mesh);
473 for (auto array_index = 0; array_index < face_data->GetNumberOfArrays(); ++array_index) {
474 auto* face_array = face_data->GetArray(array_index);
475 if (String name = face_array->GetName(); name.substring(0, 6) == "GROUP_")
476 _createGroup(face_array, name.substring(6), mesh, mesh->faceFamily(), vtk_to_Arcane_lids);
477 else
478 _createVariable(face_array, name, mesh, mesh->faceFamily(), arcane_to_vtk_lids);
479 if (m_print_info_level.print_debug_info) {
480 debug(Trace::High) << "Reading property " << face_array->GetName();
481 for (auto tuple_index = 0; tuple_index < face_array->GetNumberOfTuples(); ++tuple_index) {
482 for (auto component_index = 0; component_index < face_array->GetNumberOfComponents(); ++component_index) {
483 debug(Trace::High) << face_array->GetName() << "[" << tuple_index << "][" << component_index << "] = " << face_array->GetComponent(tuple_index, component_index);
484 }
485 }
486 }
487 }
488 }
489}
490
491/*---------------------------------------------------------------------------*/
492/*---------------------------------------------------------------------------*/
493
494void VtkPolyhedralMeshIOService::
495_createGroup(vtkDataArray* group_items, const String& group_name, IPrimaryMesh* mesh, IItemFamily* item_family, Int32ConstSpan vtkToArcaneLid) const
496{
497 ARCANE_CHECK_POINTER(group_items);
499 ARCANE_CHECK_POINTER(item_family);
500 if (group_items->GetNumberOfComponents() != 1)
501 fatal() << String::format("Cannot create item group {0}. Group information in data property must be a scalar", group_name);
502 debug() << "Create group " << group_name;
503 Int32UniqueArray arcane_lids;
504 arcane_lids.reserve((int)group_items->GetNumberOfValues());
505 using GroupDispatcher = vtkArrayDispatch::DispatchByValueType<vtkArrayDispatch::Integrals>;
506 auto group_creator = [&arcane_lids, &vtkToArcaneLid](auto* array) {
507 vtkIdType numTuples = array->GetNumberOfTuples();
508 vtkDataArrayAccessor<std::remove_pointer_t<decltype(array)>> array_accessor{ array };
509 auto local_id = 0;
510 for (vtkIdType tupleIdx = 0; tupleIdx < numTuples; ++tupleIdx) {
511 auto value = array_accessor.Get(tupleIdx, 0);
512 if (value)
513 arcane_lids.push_back(vtkToArcaneLid[local_id]);
514 ++local_id;
515 }
516 };
517 if (!GroupDispatcher::Execute(group_items, group_creator))
518 ARCANE_FATAL("Cannot create item group {0}. Group information in data property must be an integral type", group_name);
519 debug() << " local ids for item group " << group_name << " " << arcane_lids; // to remove
520 item_family->createGroup(group_name, arcane_lids);
521}
522
523/*---------------------------------------------------------------------------*/
524/*---------------------------------------------------------------------------*/
525
526template <typename vtkType>
528{
529 using type = vtkType;
530};
531
532template <>
534{
535 using type = Real;
536};
537
538template <>
540{
541 using type = Int64;
542};
543
544template <typename T> using to_arcane_type_t = typename ToArcaneType<T>::type;
545/*---------------------------------------------------------------------------*/
546
547void VtkPolyhedralMeshIOService::
549{
553 if (item_values->GetNumberOfTuples() != item_family->nbItem())
554 ARCANE_FATAL("Cannot create variable {0}, {1} values are given for {2} items in {3} family",
555 variable_name, item_values->GetNumberOfTuples(), item_family->nbItem(), item_family->name());
556 debug(Trace::High) << "Create mesh variable " << variable_name;
557 auto variable_creator = [mesh, variable_name, item_family, arcane_to_vtk_lids, this](auto* values) {
559 using ValueType = typename std::remove_pointer_t<decltype(values)>::ValueType;
560 auto* var = new ItemVariableScalarRefT<to_arcane_type_t<ValueType>>{ vbi, item_family->itemKind() };
561 m_read_variables.add(var);
562 vtkDataArrayAccessor<std::remove_pointer_t<decltype(values)>> values_accessor{ values };
563 ENUMERATE_ITEM (item, item_family->allItems()) {
564 (*var)[item] = (to_arcane_type_t<ValueType>)values_accessor.Get(arcane_to_vtk_lids[item.localId()], 0);
565 }
566 };
567 auto array_variable_creator = [mesh, variable_name, item_family, arcane_to_vtk_lids, this](auto* values) {
568 VariableBuildInfo vbi{ mesh, variable_name };
569 using ValueType = typename std::remove_pointer_t<decltype(values)>::ValueType;
570 auto* var = new ItemVariableArrayRefT<to_arcane_type_t<ValueType>>{ vbi, item_family->itemKind() };
571 m_read_variables.add(var);
572 vtkDataArrayAccessor<std::remove_pointer_t<decltype(values)>> values_accessor{ values };
573 var->resize(values->GetNumberOfComponents());
574 ENUMERATE_ITEM (item, item_family->allItems()) {
575 auto index = 0;
576 for (auto& var_value : (*var)[item]) {
577 var_value = (to_arcane_type_t<ValueType>)values_accessor.Get(arcane_to_vtk_lids[item.localId()], index++);
578 }
579 }
580 };
581 // Restrict to int and real values
582 using ValueTypes = vtkTypeList_Create_6(double, float, int, long, long long, short);
583 using ArrayDispatcher = vtkArrayDispatch::DispatchByValueType<ValueTypes>;
584 // Create ScalarVariable
585 bool is_variable_created = false;
586 if (item_values->GetNumberOfComponents() == 1) {
587 is_variable_created = ArrayDispatcher::Execute(item_values, variable_creator);
588 }
589 // Create ArrayVariable
590 else { // ArrayVariable
591 is_variable_created = ArrayDispatcher::Execute(item_values, array_variable_creator);
592 }
593 if (!is_variable_created)
594 ARCANE_FATAL("Cannot create variable {0}, it's data type is not supported. Only real and integral types are supported", variable_name);
595}
596
597/*---------------------------------------------------------------------------*/
598/*---------------------------------------------------------------------------*/
599
600void VtkPolyhedralMeshIOService::
601_computeFaceVtkArcaneLidConversion(Int32Span face_vtk_to_arcane_lids, Int32Span arcane_to_vtk_lids, VtkPolyhedralMeshIOService::VtkReader& reader, IPrimaryMesh* mesh) const
602{
603 auto face_nodes_unique_ids = reader.faceNodes();
604 auto face_nb_nodes = reader.faceNbNodes();
605 auto current_face_index = 0;
606 auto current_face_index_in_face_nodes = 0;
607 Int32UniqueArray face_nodes_local_ids(face_nodes_unique_ids.size());
608 mesh->nodeFamily()->itemsUniqueIdToLocalId(face_nodes_local_ids, face_nodes_unique_ids);
609 for (auto current_face_nb_node : face_nb_nodes) {
610 auto current_face_nodes = face_nodes_local_ids.subConstView(current_face_index_in_face_nodes, current_face_nb_node);
611 current_face_index_in_face_nodes += current_face_nb_node;
612 Node face_first_node{ mesh->nodeFamily()->view()[current_face_nodes[0]] };
613 face_vtk_to_arcane_lids[current_face_index] = mesh_utils::getFaceFromNodesLocal(face_first_node, current_face_nodes).localId();
614 ++current_face_index;
615 }
616 auto vtk_lid = 0;
617 for (auto arcane_lid : face_vtk_to_arcane_lids) {
618 arcane_to_vtk_lids[arcane_lid] = vtk_lid;
619 ++vtk_lid;
620 }
621}
622
623/*---------------------------------------------------------------------------*/
624/*---------------------------------------------------------------------------*/
625
626VtkPolyhedralMeshIOService::VtkReader::
627VtkReader(const String& filename, VtkPolyhedralTools::PrintInfoLevel print_info_level)
628: m_filename{ filename }
629, m_print_info_level{ print_info_level }
630{
631 if (filename.empty()) {
632 m_read_status.failure = true;
633 m_read_status.failure_message = "filename for polyhedral vtk mesh is empty.";
634 return;
635 }
636 if (filename.endsWith("vtk"))
637 _readPlainTextVtkGrid(filename);
638 else if (filename.endsWith("vtu"))
639 _readXlmVtkGrid(filename);
640 else {
641 m_read_status.failure = true;
642 m_read_status.failure_message = String::format("Unsupported vtk extension for file {0}. Supported vtk extension for Polyhedral meshes are {1}",
643 filename,VtkReader::supportedVtkExtensions());
644 }
645 // Check vtk grid exists and not empty
646 if (!m_vtk_grid) {
647 m_read_status.failure = true;
648 m_read_status.failure_message = String::format("Cannot read vtk polyhedral file {0}. Vtk grid was not created.", filename);
649 return;
650 }
651 if (m_vtk_grid->GetNumberOfCells() == 0) {
652 m_read_status.failure = true;
653 m_read_status.failure_message = String::format("Cannot read vtk polyhedral file {0}. No cells were found.", filename);
654 return;
655 }
656 if (!m_vtk_grid->GetFaces()) {
657 m_read_status.failure = true;
658 m_read_status.failure_message = String::format("The given mesh vtk file {0} is not a polyhedral mesh, cannot read it", filename);
659 return;
660 }
661
662 m_cell_data = m_vtk_grid->GetCellData();
663 m_point_data = m_vtk_grid->GetPointData();
664
665 // Read face info (for variables and groups) if present
666 String faces_filename = m_filename + "faces.vtk";
667 std::ifstream ifile(faces_filename.localstr());
668 if (ifile) {
669 _readPlainTextVtkFaceGrid(faces_filename);
670 }
671 else{
672 faces_filename = m_filename + "faces.vtu";
673 ifile = std::ifstream{ faces_filename.localstr() };
674 if (ifile) {
675 _readXmlVtkFaceGrid(faces_filename);
676 }
677 }
678
679 StringUniqueArray faces_filename_and_extension;
680 faces_filename.split(faces_filename_and_extension,'.');
681
682 if (!ifile){
683 m_read_status.info_message = String::format("Information no face mesh given {0}{1} (.vtk or .vtu) to define face variables or groups on faces.",
684 faces_filename_and_extension[0],
685 faces_filename_and_extension[1]);
686 }
687 else {
688 if (m_vtk_face_grid) { // Check face vtk grid exists and not empty
689 if (m_vtk_face_grid->GetNumberOfCells() == 0) {
690 m_read_status.failure = true;
691 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);
692 }
693 else {
694 m_face_data = m_vtk_face_grid->GetCellData();
695 }
696 }
697 else{
698 m_read_status.failure = true;
699 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).",
700 faces_filename_and_extension[0],
701 faces_filename_and_extension[1]);
702 }
703 }
704
705 if (m_print_info_level.print_mesh_info)
706 _printMeshInfos();
707}
708
709/*---------------------------------------------------------------------------*/
710/*---------------------------------------------------------------------------*/
711
712Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
713cellUids()
714{
715 if (m_cell_uids.empty()) {
716 _checkVtkGrid();
717 m_cell_uids.reserve(m_vtk_grid->GetNumberOfCells());
718 m_cell_nb_nodes.reserve(m_vtk_grid->GetNumberOfCells());
719 m_cell_node_uids.reserve(10 * m_vtk_grid->GetNumberOfCells()); // take a mean of 10 nodes per cell
720 auto* cell_iter = m_vtk_grid->NewCellIterator();
721 cell_iter->InitTraversal();
722 while (!cell_iter->IsDoneWithTraversal()) {
723 m_cell_uids.push_back(cell_iter->GetCellId());
724 m_cell_nb_nodes.push_back(Integer(cell_iter->GetNumberOfPoints()));
725 ArrayView<vtkIdType> cell_nodes{ Integer(cell_iter->GetNumberOfPoints()), cell_iter->GetPointIds()->GetPointer(0) };
726 std::for_each(cell_nodes.begin(), cell_nodes.end(), [this](auto uid) { this->m_cell_node_uids.push_back(uid); });
727 cell_iter->GoToNextCell();
728 }
729 }
730 return m_cell_uids;
731}
732
733/*---------------------------------------------------------------------------*/
734/*---------------------------------------------------------------------------*/
735
736Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
737nodeUids()
738{
739 if (m_node_uids.empty()) {
740 _checkVtkGrid();
741 auto nb_nodes = m_vtk_grid->GetNumberOfPoints();
742 m_node_uids.resize(nb_nodes);
743 m_node_nb_cells.resize(nb_nodes);
744 m_node_cell_uids.reserve(8 * nb_nodes);
745 m_node_uid_to_index.resize(nb_nodes);
746 for (int node_index = 0; node_index < nb_nodes; ++node_index) {
747 Int64 node_uid = node_index;
748 m_node_uids[node_index] = node_uid;
749 auto cell_nodes = vtkIdList::New();
750 m_vtk_grid->GetPointCells(node_index, cell_nodes);
751 Int64Span cell_nodes_view((Int64*)cell_nodes->GetPointer(0), cell_nodes->GetNumberOfIds());
752 m_node_cell_uids.addRange(cell_nodes_view);
753 m_node_nb_cells[node_index] = (Int32)cell_nodes->GetNumberOfIds();
754 // uid and index might differ in the future (ex parallel read).
755 // This structure is created to be used in faceUids. Should be a map if node_uid is no longer the node index
756 m_node_uid_to_index[node_uid] = node_index;
757 }
758 }
759 return m_node_uids;
760}
761
762/*---------------------------------------------------------------------------*/
763/*---------------------------------------------------------------------------*/
764
765Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
766faceUids()
767{
768 // Needs nodeUids to be called first (no work done if already called)
769 nodeUids();
770
771 if (!m_face_uids.empty())
772 return m_face_uids;
773 _checkVtkGrid();
774 auto* cell_iter = m_vtk_grid->NewCellIterator();
775 cell_iter->InitTraversal();
776 vtkIdType nb_face_estimation = 0;
777 while (!cell_iter->IsDoneWithTraversal()) {
778 vtkIdType cell_nb_faces = 0;
779 vtkIdType* points{ nullptr };
780 m_vtk_grid->GetFaceStream(cell_iter->GetCellId(), cell_nb_faces, points);
781 nb_face_estimation += cell_nb_faces;
782 cell_iter->GoToNextCell();
783 }
784 m_face_uids.reserve(nb_face_estimation);
785 auto const* faces = m_vtk_grid->GetFaces();
786 // This array contains the face info per cells (cf. vtk file)
787 // 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
788
789 if (!faces) {
790 ARCANE_FATAL("Mesh {0} is not polyhedral: faces are not defined", m_filename);
791 }
792 Int64 face_uid = 0;
793 auto face_info_size = faces->GetNumberOfValues();
794 m_face_node_uids.reserve(face_info_size);
795 m_face_nb_nodes.reserve(nb_face_estimation);
796 m_face_cell_uids.reserve(2 * nb_face_estimation);
797 m_face_nb_cells.reserve(nb_face_estimation);
798 m_cell_face_uids.reserve(8 * m_cell_uids.size()); // take a mean of 8 faces per cell
799 m_cell_nb_faces.resize(m_cell_uids.size(), 0);
800 m_cell_face_indexes.resize(m_cell_uids.size(), -1);
801 m_face_uid_indexes.resize(2 * nb_face_estimation, -1);
802 Int64UniqueArray current_face_nodes, sorted_current_face_nodes;
803 current_face_nodes.reserve(10);
804 sorted_current_face_nodes.reserve(10);
805 UniqueArray<UniqueArray<Int64>> node_faces(m_node_uids.size());
806 UniqueArray<Int32> face_offsets;
807 face_offsets.reserve(nb_face_estimation);
808 face_offsets.push_back(0);
809 FaceUidToIndexMap face_uid_to_index;
810 face_uid_to_index.reserve(nb_face_estimation);
811 auto cell_index = 0;
812 auto cell_face_index = 0;
813 auto global_face_index = 0;
814 auto face_uid_index = 0;
815 for (int face_info_index = 0; face_info_index < face_info_size; cell_index++) { // face data are given by cell
816 auto current_cell_nb_faces = Int32(faces->GetValue(face_info_index++));
817 m_cell_face_indexes[m_cell_uids[cell_index]] = cell_face_index;
818 for (auto face_index = 0; face_index < current_cell_nb_faces; ++face_index, ++global_face_index) {
819 auto current_face_nb_nodes = Int32(faces->GetValue(face_info_index++));
820 m_cell_nb_faces[m_cell_uids[cell_index]] += 1;
821 for (int node_index = 0; node_index < current_face_nb_nodes; ++node_index) {
822 current_face_nodes.push_back(faces->GetValue(face_info_index++));
823 }
824 sorted_current_face_nodes.resize(current_face_nodes.size());
825 auto is_front_cell = mesh_utils::reorderNodesOfFace(current_face_nodes, sorted_current_face_nodes);
826 auto [is_face_found, existing_face_index] = _findFace(sorted_current_face_nodes, node_faces,
827 m_node_uid_to_index, m_face_nb_nodes,
828 face_uid_to_index, face_offsets, m_face_node_uids);
829 if (!is_face_found) {
830 for (auto node_uid : current_face_nodes) {
831 node_faces[m_node_uid_to_index[node_uid]].push_back(face_uid);
832 }
833 m_cell_face_uids.push_back(face_uid);
834 m_face_uids.push_back(face_uid); // todo parallel
835 m_face_nb_nodes.push_back(current_face_nb_nodes);
836 m_face_node_uids.addRange(sorted_current_face_nodes);
837 m_face_nb_cells.push_back(1);
838 m_face_uid_indexes[global_face_index] = face_uid_index;
839 face_uid_to_index.push_back(face_uid_index);
840 auto previous_offset = face_offsets.back();
841 face_offsets.push_back(previous_offset + sorted_current_face_nodes.size());
842 ++face_uid;
843 ++face_uid_index;
844 if (is_front_cell) {
845 m_face_cell_uids.push_back(NULL_ITEM_UNIQUE_ID);
846 m_face_cell_uids.push_back(m_cell_uids[cell_index]);
847 }
848 else {
849 m_face_cell_uids.push_back(m_cell_uids[cell_index]);
850 m_face_cell_uids.push_back(NULL_ITEM_UNIQUE_ID);
851 }
852 }
853 else {
854 m_cell_face_uids.push_back(m_face_uids[existing_face_index]);
855 m_face_nb_cells[existing_face_index] += 1;
856 m_face_uid_indexes[global_face_index] = existing_face_index;
857 // add cell to face cell connectivity
858 if (is_front_cell) {
859 if (m_face_cell_uids[2 * existing_face_index + 1] != NULL_ITEM_UNIQUE_ID) {
860 ARCANE_FATAL("Problem in face orientation, face uid {0}, nodes {1}, same orientation in cell {2} and {3}. Change mesh file.",
861 m_face_uids[existing_face_index],
862 current_face_nodes,
863 m_face_cell_uids[2 * existing_face_index + 1],
864 m_cell_uids[cell_index]);
865 }
866 m_face_cell_uids[2 * existing_face_index + 1] = m_cell_uids[cell_index];
867 }
868 else {
869 if (m_face_cell_uids[2 * existing_face_index] != NULL_ITEM_UNIQUE_ID) {
870 ARCANE_FATAL("Problem in face orientation, face uid {0}, nodes {1}, same orientation in cell {2} and {3}. Change mesh file.",
871 m_face_uids[existing_face_index],
872 current_face_nodes,
873 m_face_cell_uids[2 * existing_face_index],
874 m_cell_uids[cell_index]);
875 }
876 m_face_cell_uids[2 * existing_face_index] = m_cell_uids[cell_index];
877 }
878 }
879 current_face_nodes.clear();
880 sorted_current_face_nodes.clear();
881 }
882 cell_face_index += m_cell_nb_faces[m_cell_uids[cell_index]];
883 }
884 // fill node_face_uids and node_nb_faces from node_faces (array form [nb_nodes][nb_connected_faces])
885 m_node_nb_faces.resize(m_node_uids.size(), 0);
886 _flattenConnectivity(node_faces.constSpan(), m_node_nb_faces, m_node_face_uids);
887
888 if (m_print_info_level.print_debug_info) {
889 std::cout << "================FACE NODES ==============" << std::endl;
890 std::copy(m_face_node_uids.begin(), m_face_node_uids.end(), std::ostream_iterator<Int64>(std::cout, " "));
891 std::cout << std::endl;
892 std::copy(m_face_nb_nodes.begin(), m_face_nb_nodes.end(), std::ostream_iterator<Int64>(std::cout, " "));
893 std::cout << std::endl;
894 std::copy(m_cell_face_indexes.begin(), m_cell_face_indexes.end(), std::ostream_iterator<Int64>(std::cout, " "));
895 std::cout << std::endl;
896 }
897
898 return m_face_uids;
899}
900
901/*---------------------------------------------------------------------------*/
902/*---------------------------------------------------------------------------*/
903
904Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
905edgeUids()
906{
907 if (!m_edge_uids.empty())
908 return m_edge_uids;
909
910 nodeUids(); // Needed to be called first. No work done if already called
911
912 _checkVtkGrid();
913 m_edge_uids.reserve(2 * m_vtk_grid->GetNumberOfPoints());
914 auto const* faces = m_vtk_grid->GetFaces();
915 // This array contains the face info per cells (cf. vtk file)
916 // 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
917
918 if (!faces) {
919 ARCANE_FATAL("Mesh {0} is not polyhedral: faces are not defined", m_filename);
920 }
921 Int64 edge_uid = 0;
922 auto nb_edge_estimation = 2 * m_edge_uids.capacity();
923 m_edge_node_uids.reserve(nb_edge_estimation);
924 auto face_info_size = faces->GetNumberOfValues();
925 auto cell_index = 0;
926 auto global_face_index = 0;
927 auto new_edge_index = 0;
928 UniqueArray<std::set<Int64>> edge_cells;
929 UniqueArray<Int64UniqueArray> edge_faces;
930 edge_cells.reserve(m_edge_uids.capacity());
931 edge_faces.reserve(m_edge_uids.capacity());
932 m_cell_nb_edges.resize(m_cell_uids.size(), 0);
933 m_cell_edge_uids.reserve(20 * m_cell_uids.size()); // choose a value of 20 edge per cell
934 UniqueArray<std::set<Int64>> face_edges;
935 face_edges.resize(m_face_uids.size());
936 UniqueArray<std::set<Int64>> cell_edges;
937 cell_edges.resize(m_cell_uids.size());
938 UniqueArray<Int64UniqueArray> node_edges;
939 node_edges.resize(m_node_uids.size());
940 EdgeUidToIndexMap edge_uid_to_index;
941 edge_uid_to_index.reserve(nb_edge_estimation);
942 Int32UniqueArray edge_offsets;
943 edge_offsets.reserve(nb_edge_estimation);
944 edge_offsets.push_back(0);
945 m_edge_nb_nodes.reserve(nb_edge_estimation);
946 for (int face_info_index = 0; face_info_index < face_info_size; ++cell_index) {
947 auto current_cell_nb_faces = Int32(faces->GetValue(face_info_index++));
948 for (auto face_index = 0; face_index < current_cell_nb_faces; ++face_index, ++global_face_index) {
949 auto current_face_nb_nodes = Int32(faces->GetValue(face_info_index++));
950 auto first_face_node_uid = Int32(faces->GetValue(face_info_index));
951 UniqueArray<Int64> current_edge(2), sorted_edge(2);
952 for (int node_index = 0; node_index < current_face_nb_nodes - 1; ++node_index) {
953 current_edge = UniqueArray<Int64>{ faces->GetValue(face_info_index++), faces->GetValue(face_info_index) };
954 mesh_utils::reorderNodesOfFace(current_edge, sorted_edge); // works for edges
955 auto [is_edge_found, existing_edge_index] = _findFace(sorted_edge, node_edges,
956 m_node_uid_to_index,
957 m_edge_nb_nodes,
958 edge_uid_to_index, edge_offsets, m_edge_node_uids); // works for edges
959 if (!is_edge_found) {
960 m_cell_nb_edges[cell_index] += 1;
961 face_edges[m_face_uid_indexes[global_face_index]].insert(edge_uid);
962 cell_edges[cell_index].insert(edge_uid);
963 for (auto node : current_edge) {
964 node_edges[node].push_back(edge_uid);
965 }
966 edge_cells.push_back(std::set{ m_cell_uids[cell_index] });
967 edge_faces.push_back(Int64UniqueArray{ m_cell_face_uids[m_cell_face_indexes[cell_index] + face_index] });
968 m_edge_uids.push_back(edge_uid++); // todo parallel
969 m_edge_node_uids.addRange(sorted_edge);
970 edge_uid_to_index.push_back(new_edge_index);
971 auto current_offset = edge_offsets.back();
972 edge_offsets.push_back(current_offset + 2);
973 m_edge_nb_nodes.push_back(2);
974 ++new_edge_index;
975 }
976 else {
977 edge_cells[existing_edge_index].insert(m_cell_uids[cell_index]);
978 edge_faces[existing_edge_index].push_back(m_cell_face_uids[m_cell_face_indexes[cell_index] + face_index]);
979 face_edges[m_face_uid_indexes[global_face_index]].insert(m_edge_uids[existing_edge_index]);
980 cell_edges[cell_index].insert(m_edge_uids[existing_edge_index]);
981 }
982 }
983 current_edge = UniqueArray<Int64>{ faces->GetValue(face_info_index++), first_face_node_uid };
984 mesh_utils::reorderNodesOfFace(current_edge, sorted_edge); // works for edges
985 auto [is_edge_found, existing_edge_index] = _findFace(sorted_edge, node_edges,
986 m_node_uid_to_index,
987 m_edge_nb_nodes,
988 edge_uid_to_index, edge_offsets, m_edge_node_uids); // works for edges
989 if (!is_edge_found) {
990 m_cell_nb_edges[cell_index] += 1;
991 edge_cells.push_back(std::set{ m_cell_uids[cell_index] });
992 edge_faces.push_back(Int64UniqueArray{ m_cell_face_uids[m_cell_face_indexes[cell_index] + face_index] });
993 face_edges[m_face_uid_indexes[global_face_index]].insert(edge_uid);
994 cell_edges[cell_index].insert(edge_uid);
995 for (auto node : current_edge) {
996 node_edges[node].push_back(edge_uid);
997 }
998 m_edge_uids.push_back(edge_uid++); // todo parallel
999 m_edge_node_uids.addRange(sorted_edge);
1000 edge_uid_to_index.push_back(new_edge_index);
1001 auto current_offset = edge_offsets.back();
1002 edge_offsets.push_back(current_offset + 2);
1003 m_edge_nb_nodes.push_back(2);
1004 ++new_edge_index;
1005 }
1006 else {
1007 edge_cells[existing_edge_index].insert(m_cell_uids[cell_index]);
1008 edge_faces[existing_edge_index].push_back(m_cell_face_uids[m_cell_face_indexes[cell_index] + face_index]);
1009 face_edges[m_face_uid_indexes[global_face_index]].insert(m_edge_uids[existing_edge_index]);
1010 cell_edges[cell_index].insert(m_edge_uids[existing_edge_index]);
1011 }
1012 }
1013 }
1014 // fill edge_cell_uids and edge_nb_cells from edge_cells (array form [nb_edges][nb_connected_cells])
1015 m_edge_nb_cells.resize(m_edge_uids.size(), 0);
1016 _flattenConnectivity(edge_cells.constSpan(), m_edge_nb_cells, m_edge_cell_uids);
1017
1018 // fill edge faces uids
1019 m_edge_nb_faces.resize(m_edge_uids.size(), 0);
1020 _flattenConnectivity(edge_faces.constSpan(), m_edge_nb_faces, m_edge_face_uids);
1021
1022 // fill face edge uids
1023 m_face_nb_edges.resize(m_face_uids.size(), 0);
1024 _flattenConnectivity(face_edges.constSpan(), m_face_nb_edges, m_face_edge_uids);
1025
1026 // fill cell edge uids
1027 m_cell_nb_edges.resize(m_cell_uids.size(), 0);
1028 _flattenConnectivity(cell_edges, m_cell_nb_edges, m_cell_edge_uids);
1029
1030 // fill node edge uids
1031 m_node_nb_edges.resize(m_node_uids.size(), 0);
1032 _flattenConnectivity(node_edges, m_node_nb_edges, m_node_edge_uids);
1033
1034 if (m_print_info_level.print_debug_info) {
1035 std::cout << "================EDGE NODES ==============" << std::endl;
1036 std::copy(m_edge_node_uids.begin(), m_edge_node_uids.end(), std::ostream_iterator<Int64>(std::cout, " "));
1037 std::cout << std::endl;
1038 std::cout << "================FACE EDGES ==============" << std::endl;
1039 std::copy(m_face_nb_edges.begin(), m_face_nb_edges.end(), std::ostream_iterator<Int32>(std::cout, " "));
1040 std::cout << std::endl;
1041 std::copy(m_face_edge_uids.begin(), m_face_edge_uids.end(), std::ostream_iterator<Int64>(std::cout, " "));
1042 std::cout << std::endl;
1043 std::cout << "================CELL EDGES ==============" << std::endl;
1044 std::copy(m_cell_nb_edges.begin(), m_cell_nb_edges.end(), std::ostream_iterator<Int32>(std::cout, " "));
1045 std::cout << std::endl;
1046 std::copy(m_cell_edge_uids.begin(), m_cell_edge_uids.end(), std::ostream_iterator<Int64>(std::cout, " "));
1047 std::cout << std::endl;
1048 }
1049 return m_edge_uids;
1050}
1051
1052/*---------------------------------------------------------------------------*/
1053/*---------------------------------------------------------------------------*/
1054
1055std::pair<bool, Int32> VtkPolyhedralMeshIOService::VtkReader::
1056_findFace(const Int64UniqueArray& sorted_face_nodes,
1057 const UniqueArray<Int64UniqueArray>& node_face_uids,
1058 const NodeUidToIndexMap& node_uid_to_index,
1059 const Int32UniqueArray& face_nb_nodes,
1060 const FaceUidToIndexMap& face_uid_to_index,
1061 const UniqueArray<Int32>& face_offsets,
1062 const Int64UniqueArray& face_node_uids)
1063{
1064 auto first_node_uid = sorted_face_nodes[0];
1065 auto first_node_index = node_uid_to_index[first_node_uid];
1066 // If the face already exists it has already been registered in node_face connectivity
1067 for (auto face_uid : node_face_uids[first_node_index]) {
1068 auto face_index = face_uid_to_index[face_uid];
1069 auto face_offset = face_offsets[face_index];
1070 auto face_nb_node = face_nb_nodes[face_index];
1071 if (face_nb_node == sorted_face_nodes.size()) {
1072 bool is_same_face = true;
1073 for (auto index = 0; index < face_nb_node; ++index) {
1074 if (sorted_face_nodes[index] != face_node_uids[face_offset + index]) {
1075 is_same_face = false;
1076 }
1077 }
1078 if (is_same_face)
1079 return { true, face_index };
1080 }
1081 }
1082 return { false, -1 };
1083}
1084
1085/*---------------------------------------------------------------------------*/
1086/*---------------------------------------------------------------------------*/
1087
1088Integer VtkPolyhedralMeshIOService::VtkReader::
1089nbNodes()
1090{
1091 if (m_node_uids.empty())
1092 nodeUids();
1093 return m_node_uids.size();
1094}
1095
1096/*---------------------------------------------------------------------------*/
1097/*---------------------------------------------------------------------------*/
1098
1099Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1100cellNodes()
1101{
1102 if (m_cell_node_uids.empty())
1103 cellUids();
1104 return m_cell_node_uids;
1105}
1106
1107/*---------------------------------------------------------------------------*/
1108/*---------------------------------------------------------------------------*/
1109
1110Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1111cellNbNodes()
1112{
1113 if (m_cell_nb_nodes.empty())
1114 cellUids();
1115 return m_cell_nb_nodes;
1116}
1117
1118/*---------------------------------------------------------------------------*/
1119/*---------------------------------------------------------------------------*/
1120
1121Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1122faceNodes()
1123{
1124 if (m_face_node_uids.empty())
1125 faceUids();
1126 return m_face_node_uids;
1127}
1128
1129/*---------------------------------------------------------------------------*/
1130/*---------------------------------------------------------------------------*/
1131
1132Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1133faceNbNodes()
1134{
1135 if (m_face_nb_nodes.empty())
1136 faceUids();
1137 return m_face_nb_nodes;
1138}
1139
1140/*---------------------------------------------------------------------------*/
1141/*---------------------------------------------------------------------------*/
1142
1143Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1144edgeNbNodes()
1145{
1146 if (m_edge_node_uids.empty())
1147 edgeUids();
1148 return m_edge_nb_nodes;
1149}
1150
1151/*---------------------------------------------------------------------------*/
1152/*---------------------------------------------------------------------------*/
1153
1154Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1155edgeNodes()
1156{
1157 if (m_edge_node_uids.empty())
1158 edgeUids();
1159 return m_edge_node_uids;
1160}
1161
1162/*---------------------------------------------------------------------------*/
1163/*---------------------------------------------------------------------------*/
1164
1165Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1166faceCells()
1167{
1168 if (m_face_cell_uids.empty())
1169 faceUids();
1170 // debug
1171 if (m_print_info_level.print_debug_info) {
1172 std::cout << "=================FACE CELLS================="
1173 << "\n";
1174 std::copy(m_face_cell_uids.begin(), m_face_cell_uids.end(), std::ostream_iterator<Int64>(std::cout, " "));
1175 std::cout << "\n";
1176 std::cout << "=================END FACE CELLS================="
1177 << "\n";
1178 }
1179 return m_face_cell_uids;
1180}
1181
1182/*---------------------------------------------------------------------------*/
1183/*---------------------------------------------------------------------------*/
1184
1185Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1186faceNbCells()
1187{
1188 if (m_face_nb_cells.empty())
1189 faceUids();
1190 return m_face_nb_cells;
1191}
1192
1193/*---------------------------------------------------------------------------*/
1194/*---------------------------------------------------------------------------*/
1195
1196Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1197edgeNbCells()
1198{
1199 if (m_edge_nb_cells.empty())
1200 edgeUids();
1201 return m_edge_nb_cells;
1202}
1203
1204/*---------------------------------------------------------------------------*/
1205/*---------------------------------------------------------------------------*/
1206
1207Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1208edgeCells()
1209{
1210 if (m_edge_cell_uids.empty())
1211 edgeUids();
1212 return m_edge_cell_uids;
1213}
1214
1215/*---------------------------------------------------------------------------*/
1216/*---------------------------------------------------------------------------*/
1217
1218Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1219cellNbFaces()
1220{
1221 if (m_cell_nb_faces.empty())
1222 faceUids();
1223 return m_cell_nb_faces;
1224}
1225
1226/*---------------------------------------------------------------------------*/
1227/*---------------------------------------------------------------------------*/
1228
1229Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1230cellFaces()
1231{
1232 if (m_cell_face_uids.empty())
1233 faceUids();
1234 return m_cell_face_uids;
1235}
1236
1237/*---------------------------------------------------------------------------*/
1238/*---------------------------------------------------------------------------*/
1239
1240Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1241edgeNbFaces()
1242{
1243 if (m_edge_nb_faces.empty())
1244 edgeUids();
1245 return m_edge_nb_faces;
1246}
1247
1248/*---------------------------------------------------------------------------*/
1249/*---------------------------------------------------------------------------*/
1250
1251Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1252edgeFaces()
1253{
1254 if (m_edge_face_uids.empty())
1255 edgeUids();
1256 return m_edge_face_uids;
1257}
1258
1259/*---------------------------------------------------------------------------*/
1260/*---------------------------------------------------------------------------*/
1261
1262Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1263cellNbEdges()
1264{
1265 if (m_cell_nb_edges.empty())
1266 edgeUids();
1267 return m_cell_nb_edges;
1268}
1269
1270/*---------------------------------------------------------------------------*/
1271/*---------------------------------------------------------------------------*/
1272
1273Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1274cellEdges()
1275{
1276 if (m_cell_edge_uids.empty())
1277 edgeUids();
1278 return m_cell_edge_uids;
1279}
1280
1281/*---------------------------------------------------------------------------*/
1282/*---------------------------------------------------------------------------*/
1283
1284Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1285faceNbEdges()
1286{
1287 if (m_face_nb_edges.empty())
1288 edgeUids();
1289 return m_face_nb_edges;
1290}
1291
1292/*---------------------------------------------------------------------------*/
1293/*---------------------------------------------------------------------------*/
1294
1295Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1296faceEdges()
1297{
1298 if (m_face_edge_uids.empty())
1299 edgeUids();
1300 return m_face_edge_uids;
1301}
1302
1303/*---------------------------------------------------------------------------*/
1304/*---------------------------------------------------------------------------*/
1305
1306template <typename Connectivity2DArray>
1307void VtkPolyhedralMeshIOService::VtkReader::
1308_flattenConnectivity(Connectivity2DArray connected_item_2darray,
1309 Int32Span nb_connected_item_per_source_item,
1310 Int64UniqueArray& connected_item_array)
1311{
1312 // fill nb_connected_item_per_source_items
1313 std::transform(connected_item_2darray.begin(), connected_item_2darray.end(), nb_connected_item_per_source_item.begin(), [](auto const& connected_items) {
1314 return connected_items.size();
1315 });
1316 // fill connected_item_array
1317 connected_item_array.reserve(std::accumulate(nb_connected_item_per_source_item.begin(), nb_connected_item_per_source_item.end(), 0));
1318 std::for_each(connected_item_2darray.begin(), connected_item_2darray.end(), [&connected_item_array](auto const& connected_items) {
1319 for (auto const& connected_item : connected_items) {
1320 connected_item_array.push_back(connected_item);
1321 }
1322 });
1323}
1324
1325/*---------------------------------------------------------------------------*/
1326/*---------------------------------------------------------------------------*/
1327
1328Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1329nodeNbCells()
1330{
1331 if (m_node_nb_cells.empty())
1332 nodeUids();
1333 return m_node_nb_cells;
1334}
1335
1336/*---------------------------------------------------------------------------*/
1337/*---------------------------------------------------------------------------*/
1338
1339Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1340nodeCells()
1341{
1342 if (m_node_cell_uids.empty())
1343 nodeUids();
1344 return m_node_cell_uids;
1345}
1346
1347/*---------------------------------------------------------------------------*/
1348/*---------------------------------------------------------------------------*/
1349
1350Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1351nodeNbFaces()
1352{
1353 if (m_node_nb_faces.empty())
1354 faceUids();
1355 return m_node_nb_faces;
1356}
1357
1358/*---------------------------------------------------------------------------*/
1359/*---------------------------------------------------------------------------*/
1360
1361Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1362nodeFaces()
1363{
1364 if (m_node_face_uids.empty())
1365 faceUids();
1366 return m_node_face_uids;
1367}
1368
1369/*---------------------------------------------------------------------------*/
1370/*---------------------------------------------------------------------------*/
1371
1372Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1373nodeNbEdges()
1374{
1375 if (m_node_nb_edges.empty())
1376 edgeUids();
1377 return m_node_nb_edges;
1378}
1379
1380/*---------------------------------------------------------------------------*/
1381/*---------------------------------------------------------------------------*/
1382
1383Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
1384nodeEdges()
1385{
1386 if (m_node_edge_uids.empty())
1387 edgeUids();
1388 return m_node_edge_uids;
1389}
1390
1391/*---------------------------------------------------------------------------*/
1392/*---------------------------------------------------------------------------*/
1393
1394Real3ArrayView VtkPolyhedralMeshIOService::VtkReader::
1395nodeCoords()
1396{
1397 if (m_node_coordinates.empty()) {
1398 _checkVtkGrid();
1399 auto point_coords = m_vtk_grid->GetPoints()->GetData();
1400 if (m_print_info_level.print_debug_info) {
1401 std::cout << "======= Point COORDS ====" << std::endl;
1402 std::ostringstream oss;
1403 point_coords->PrintSelf(oss, vtkIndent{ 2 });
1404 std::cout << oss.str() << std::endl;
1405 }
1406 auto nb_nodes = m_vtk_grid->GetNumberOfPoints();
1407 for (int i = 0; i < nb_nodes; ++i) {
1408 if (m_print_info_level.print_debug_info) {
1409 std::cout << "==========current point coordinates : ( ";
1410 std::cout << *(point_coords->GetTuple(i)) << " , ";
1411 std::cout << *(point_coords->GetTuple(i) + 1) << " , ";
1412 std::cout << *(point_coords->GetTuple(i) + 2) << " ) ===" << std::endl;
1413 }
1414 m_node_coordinates.add({ *(point_coords->GetTuple(i)),
1415 *(point_coords->GetTuple(i) + 1),
1416 *(point_coords->GetTuple(i) + 2) });
1417 }
1418 }
1419 return m_node_coordinates;
1420}
1421
1422/*---------------------------------------------------------------------------*/
1423/*---------------------------------------------------------------------------*/
1424
1425vtkCellData* VtkPolyhedralMeshIOService::VtkReader::
1426cellData()
1427{
1428 return m_cell_data;
1429}
1430
1431/*---------------------------------------------------------------------------*/
1432/*---------------------------------------------------------------------------*/
1433
1434vtkPointData* VtkPolyhedralMeshIOService::VtkReader::
1435pointData()
1436{
1437 return m_point_data;
1438}
1439
1440/*---------------------------------------------------------------------------*/
1441/*---------------------------------------------------------------------------*/
1442
1443void VtkPolyhedralMeshIOService::VtkReader::
1444_printMeshInfos() const
1445{
1446 _checkVtkGrid();
1447 std::cout << "-- VTK GRID READ "
1448 << " NB CELLS " << m_vtk_grid->GetNumberOfCells() << std::endl;
1449 // Parse cells
1450 auto* cell_iter = m_vtk_grid->vtkDataSet::NewCellIterator();
1451 cell_iter->InitTraversal();
1452 vtkIdType* cell_faces{ nullptr };
1453 vtkIdType nb_faces = 0;
1454 while (!cell_iter->IsDoneWithTraversal()) {
1455 std::cout << "---- visiting cell id " << cell_iter->GetCellId() << std::endl;
1456 std::cout << "---- cell number of faces " << cell_iter->GetNumberOfFaces() << std::endl;
1457 std::cout << "---- cell number of points " << cell_iter->GetNumberOfPoints() << std::endl;
1458 m_vtk_grid->GetFaceStream(cell_iter->GetCellId(), nb_faces, cell_faces);
1459 for (auto iface = 0; iface < nb_faces; ++iface) {
1460 auto face_nb_nodes = *cell_faces++;
1461 std::cout << "---- has face with " << face_nb_nodes << " nodes. Node ids : ";
1462 for (int inode = 0; inode < face_nb_nodes; ++inode) {
1463 std::cout << *cell_faces++ << " ";
1464 }
1465 std::cout << std::endl;
1466 }
1467 cell_iter->GoToNextCell();
1468 }
1469}
1470
1471/*---------------------------------------------------------------------------*/
1472/*---------------------------------------------------------------------------*/
1473
1474vtkCellData* VtkPolyhedralMeshIOService::VtkReader::faceData()
1475{
1476 return m_face_data;
1477}
1478/*---------------------------------------------------------------------------*/
1479/*---------------------------------------------------------------------------*/
1480
1481/*---------------------------------------------------------------------------*/
1482/*---------------------------------------------------------------------------*/
1483
1484void VtkPolyhedralMeshIOService::VtkReader::
1485_readPlainTextVtkGrid(const String& filename)
1486{
1487 m_vtk_grid_reader->SetFileName(filename.localstr());
1488 m_vtk_grid_reader->ReadAllScalarsOn();
1489 m_vtk_grid_reader->Update();
1490 m_vtk_grid = m_vtk_grid_reader->GetOutput();
1491}
1492
1493/*---------------------------------------------------------------------------*/
1494/*---------------------------------------------------------------------------*/
1495
1496void VtkPolyhedralMeshIOService::VtkReader::
1497_readXlmVtkGrid(const String& filename)
1498{
1499 m_vtk_xml_grid_reader->SetFileName(filename.localstr());
1500 m_vtk_xml_grid_reader->Update();
1501 m_vtk_grid = m_vtk_xml_grid_reader->GetOutput();
1502}
1503
1504/*---------------------------------------------------------------------------*/
1505/*---------------------------------------------------------------------------*/
1506
1507void VtkPolyhedralMeshIOService::VtkReader::
1508_readPlainTextVtkFaceGrid(const String& faces_filename)
1509{
1510 m_vtk_face_grid_reader->SetFileName(faces_filename.localstr());
1511 m_vtk_face_grid_reader->ReadAllScalarsOn();
1512 m_vtk_face_grid_reader->Update();
1513 m_vtk_face_grid = m_vtk_face_grid_reader->GetOutput();
1514}
1515
1516/*---------------------------------------------------------------------------*/
1517/*---------------------------------------------------------------------------*/
1518
1519void VtkPolyhedralMeshIOService::VtkReader::
1520_readXmlVtkFaceGrid(const String& faces_filename)
1521{
1522 m_vtk_xml_face_grid_reader->SetFileName(faces_filename.localstr());
1523 m_vtk_xml_face_grid_reader->Update();
1524 m_vtk_face_grid = m_vtk_xml_face_grid_reader->GetOutput();
1525}
1526
1527/*---------------------------------------------------------------------------*/
1528/*---------------------------------------------------------------------------*/
1529
1530void VtkPolyhedralMeshIOService::VtkReader::
1531_checkVtkGrid() const
1532{
1533 if (!m_vtk_grid)
1534 ARCANE_FATAL("Polyhedral vtk grid not loaded. Cannot continue.");
1535}
1536
1537/*---------------------------------------------------------------------------*/
1538/*---------------------------------------------------------------------------*/
1539
1540} // 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.
Informations nécessaires pour la lecture d'un fichier de maillage.
Interface d'une famille d'entités.
Interface d'un service de création/lecture du maillage.
virtual IMeshInitialAllocator * initialAllocator()
Allocateur initial spécifique.
Lecteur des fichiers de maillage via la bibliothèque LIMA.
Definition Lima.cc:120
Paramètres nécessaires à la construction d'un maillage.
Caractéristiques d'un maillage.
Definition MeshKind.h:59
Structure contenant les informations pour créer un service.
Paramètres nécessaires à la construction d'une variable.
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...
bool empty() const
Capacité (nombre d'éléments alloués) du vecteur.
Vue modifiable d'un tableau d'un type T.
void reserve(Int64 new_capacity)
Réserve le mémoire pour new_capacity éléments.
void push_back(ConstReferenceType val)
Ajoute l'élément val à la fin du tableau.
Vue constante d'un tableau de type T.
Interface du gestionnaire de traces.
virtual TraceMessage info()=0
Flot pour un message d'information.
Vue d'un tableau d'éléments de type T.
Definition Span.h:510
Chaîne de caractères unicode.
ITraceMng * traceMng() const
Gestionnaire de trace.
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flot pour un message de debug.
TraceMessage fatal() const
Flot pour un message d'erreur fatale.
Vecteur 1D de données avec sémantique par valeur (style STL).
#define ARCANE_REGISTER_SERVICE(aclass, a_service_property,...)
Macro pour enregistrer un service.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Span< Int32 > Int32Span
Equivalent C d'un tableau à une dimension d'entiers 32 bits.
Definition UtilsTypes.h:680
Span< Int64 > Int64Span
Equivalent C d'un tableau à une dimension d'entiers 64 bits.
Definition UtilsTypes.h:678
UniqueArray< Int64 > Int64UniqueArray
Tableau dynamique à une dimension d'entiers 64 bits.
Definition UtilsTypes.h:513
ArrayView< Real3 > Real3ArrayView
Equivalent C d'un tableau à une dimension de Real3.
Definition UtilsTypes.h:625
ConstArrayView< Int32 > Int32ConstArrayView
Equivalent C d'un tableau à une dimension d'entiers 32 bits.
Definition UtilsTypes.h:640
@ ST_CaseOption
Le service s'utilise au niveau du jeu de données.
SharedArray< Int32 > Int32SharedArray
Tableau dynamique à une dimension d'entiers 32 bits.
Definition UtilsTypes.h:547
UniqueArray< Int32 > Int32UniqueArray
Tableau dynamique à une dimension d'entiers 32 bits.
Definition UtilsTypes.h:515
@ 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.
UniqueArray< String > StringUniqueArray
Tableau dynamique à une dimension de chaînes de caractères.
Definition UtilsTypes.h:525
Span< const Int32 > Int32ConstSpan
Vue en lecture seule d'un tableau à une dimension d'entiers 32 bits.
Definition UtilsTypes.h:712
ConstArrayView< Int64 > Int64ConstArrayView
Equivalent C d'un tableau à une dimension d'entiers 64 bits.
Definition UtilsTypes.h:638
Int32 Integer
Type représentant un entier.