20#include <vtkUnstructuredGrid.h>
21#include <vtkUnstructuredGridReader.h>
23#include <vtkCellIterator.h>
24#include <vtkIdTypeArray.h>
25#include <vtkCellData.h>
26#include <vtkPointData.h>
27#include <vtkDataSetAttributes.h>
28#include <vtkArrayDispatch.h>
29#include <vtkDataArrayAccessor.h>
30#include <vtkPolyDataReader.h>
31#include <vtkPolyData.h>
32#include <vtkSmartPointer.h>
33#include <vtkCellArray.h>
34#include <vtkVersion.h>
35#if VTK_VERSION_NUMBER >= 900000000
36using vtkIdType_generic = vtkIdType
const;
38using vtkIdType_generic = vtkIdType;
42#include <arccore/base/String.h>
43#include <arccore/base/FatalErrorException.h>
44#include <vtkXMLUnstructuredGridReader.h>
45#include <vtkXMLPolyDataReader.h>
48#include "arcane/core/AbstractService.h"
49#include "arcane/core/ICaseMeshReader.h"
51#include "arcane/core/IMeshBuilder.h"
52#include "arcane/core/MeshBuildInfo.h"
53#include "arcane/core/IPrimaryMesh.h"
55#include "arcane/core/IMeshInitialAllocator.h"
56#include "arcane/core/IVariableMng.h"
58#include "arcane/utils/ITraceMng.h"
59#include "arcane/utils/UniqueArray.h"
60#include "arcane/utils/Real3.h"
61#include "arcane/mesh/CellFamily.h"
62#include "arcane/core/MeshVariableScalarRef.h"
63#include "arcane/core/MeshVariableArrayRef.h"
65#include "arcane/core/ItemAllocationInfo.h"
66#include "arcane/core/VariableBuildInfo.h"
68#include "arcane/utils/OStringStream.h"
70#include "arcane/core/IXmlDocumentHolder.h"
71#include "arcane/core/XmlNode.h"
72#include "arcane/core/internal/IVariableMngInternal.h"
73#include "arcane/core/datatype/DataTypeTraits.h"
75#include "arcane/std/VtkPolyhedralMeshIO_axl.h"
86namespace VtkPolyhedralTools
97 bool print_mesh_info =
false;
98 bool print_debug_info =
false;
105class VtkPolyhedralMeshIOService
112 , m_print_info_level(print_info_level)
122 VtkReader() =
default;
124 static String supportedVtkExtensions()
noexcept {
return "vtk,vtu"; };
174 bool readHasFailed()
const noexcept {
return m_read_status.failure; }
177 vtkCellData* cellData();
178 vtkPointData* pointData();
179 vtkCellData* faceData();
181 bool isEmpty()
const noexcept {
return m_is_empty; }
182 bool doRead()
const noexcept {
return m_do_read; }
186 bool m_is_empty =
true;
187 bool m_do_read =
false;
191 vtkNew<vtkUnstructuredGridReader> m_vtk_grid_reader;
192 vtkNew<vtkXMLUnstructuredGridReader> m_vtk_xml_grid_reader;
193 vtkNew<vtkPolyDataReader> m_vtk_face_grid_reader;
194 vtkNew<vtkXMLPolyDataReader> m_vtk_xml_face_grid_reader;
195 vtkUnstructuredGrid* m_vtk_grid =
nullptr;
196 vtkPolyData* m_vtk_face_grid =
nullptr;
212 NodeUidToIndexMap m_node_uid_to_index;
214 vtkCellData* m_cell_data =
nullptr;
215 vtkPointData* m_point_data =
nullptr;
216 vtkCellData* m_face_data =
nullptr;
217 vtkCellArray* m_poly_data =
nullptr;
219 void _printMeshInfos()
const;
222 template <
typename Connectivity2DArray>
223 static void _flattenConnectivity(Connectivity2DArray connected_item_2darray,
Int32Span nb_connected_item_per_source_item,
Int64UniqueArray& connected_item_array);
224 void _readPlainTextVtkGrid(
const String& filename);
225 void _readXlmVtkGrid(
const String& filename);
226 void _checkVtkGrid()
const;
227 void _readPlainTextVtkFaceGrid(
const String& faces_filename);
228 void _readXmlVtkFaceGrid(
const String& faces_filename);
229 void _readfaceNodesInFaceMesh();
239 bool do_read = is_parallel_read ?
mesh->parallelMng()->isMasterIO() :
true;
240 using VtkReaderPtr = std::unique_ptr<VtkReader>;
241 VtkReaderPtr reader = std::make_unique<VtkReader>();
243 reader = std::make_unique<VtkReader>(filename, m_print_info_level);
244 if (reader->readHasFailed())
245 return reader->readStatus();
247 _fillItemAllocationInfo(item_allocation_info, *reader);
248 auto polyhedral_mesh_allocator =
mesh->initialAllocator()->polyhedralMeshAllocator();
249 polyhedral_mesh_allocator->allocateItems(item_allocation_info);
250 _readVariablesAndGroups(
mesh, *reader);
251 return reader->readStatus();
261 bool is_array =
false;
271 template <
template <
class>
class VariableRootType,
template <
class>
class ArrayVariableRootType>
279class VtkPolyhedralCaseMeshReader
290 , m_read_info(read_info)
291 , m_print_info_level(print_info_level)
306 m_trace_mng->info() <<
"---Create Polyhedral mesh: " << pm->
name() <<
"---";
307 m_trace_mng->info() <<
"--Read mesh file " << m_read_info.fileName();
309 auto read_status = polyhedral_vtk_service.read(pm, m_read_info.fileName(), m_read_info.isParallelRead());
310 if (read_status.failure)
312 m_trace_mng->info() << read_status.info_message;
331 if (VtkPolyhedralMeshIOService::VtkReader::supportedVtkExtensions().contains(read_info.format()))
341 ServiceProperty(
"VtkPolyhedralCaseMeshReader",
ST_CaseOption),
347void VtkPolyhedralMeshIOService::
348_fillItemAllocationInfo(ItemAllocationInfo& item_allocation_info, VtkReader& vtk_reader)
350 auto nb_item_family = 4;
351 auto nb_connected_family = 3;
352 item_allocation_info.family_infos.resize(nb_item_family);
353 for (
auto& family_info : item_allocation_info.family_infos) {
354 family_info.connected_family_infos.resize(nb_connected_family);
357 auto& cell_family_info = item_allocation_info.family_infos[0];
358 cell_family_info.name =
"Cell";
359 cell_family_info.item_kind =
IK_Cell;
360 cell_family_info.item_uids = vtk_reader.cellUids();
361 auto& node_family_info = item_allocation_info.family_infos[1];
362 node_family_info.name =
"Node";
363 node_family_info.item_kind =
IK_Node;
364 node_family_info.item_uids = vtk_reader.nodeUids();
365 auto& face_family_info = item_allocation_info.family_infos[2];
366 face_family_info.name =
"Face";
367 face_family_info.item_kind =
IK_Face;
368 face_family_info.item_uids = vtk_reader.faceUids();
369 auto& edge_family_info = item_allocation_info.family_infos[3];
370 edge_family_info.name =
"Edge";
371 edge_family_info.item_kind =
IK_Edge;
372 edge_family_info.item_uids = vtk_reader.edgeUids();
374 auto cell_connected_family_index = 0;
375 auto& cell_connected_node_family_info = cell_family_info.connected_family_infos[cell_connected_family_index++];
376 cell_connected_node_family_info.name = node_family_info.name;
377 cell_connected_node_family_info.item_kind = node_family_info.item_kind;
378 cell_connected_node_family_info.connectivity_name =
"CellToNodes";
379 cell_connected_node_family_info.nb_connected_items_per_item = vtk_reader.cellNbNodes();
380 cell_connected_node_family_info.connected_items_uids = vtk_reader.cellNodes();
382 auto& cell_connected_face_family_info = cell_family_info.connected_family_infos[cell_connected_family_index++];
383 cell_connected_face_family_info.name = face_family_info.name;
384 cell_connected_face_family_info.item_kind = face_family_info.item_kind;
385 cell_connected_face_family_info.connectivity_name =
"CellToFaces";
386 cell_connected_face_family_info.nb_connected_items_per_item = vtk_reader.cellNbFaces();
387 cell_connected_face_family_info.connected_items_uids = vtk_reader.cellFaces();
389 auto& cell_connected_edge_family_info = cell_family_info.connected_family_infos[cell_connected_family_index++];
390 cell_connected_edge_family_info.name = edge_family_info.name;
391 cell_connected_edge_family_info.item_kind = edge_family_info.item_kind;
392 cell_connected_edge_family_info.connectivity_name =
"CellToEdges";
393 cell_connected_edge_family_info.nb_connected_items_per_item = vtk_reader.cellNbEdges();
394 cell_connected_edge_family_info.connected_items_uids = vtk_reader.cellEdges();
396 auto face_connected_family_index = 0;
397 auto& face_connected_cell_family_info = face_family_info.connected_family_infos[face_connected_family_index++];
398 face_connected_cell_family_info.name = cell_family_info.name;
399 face_connected_cell_family_info.item_kind = cell_family_info.item_kind;
400 face_connected_cell_family_info.connectivity_name =
"FaceToCells";
401 face_connected_cell_family_info.nb_connected_items_per_item = vtk_reader.faceNbCells();
402 face_connected_cell_family_info.connected_items_uids = vtk_reader.faceCells();
404 auto& face_connected_node_family_info = face_family_info.connected_family_infos[face_connected_family_index++];
405 face_connected_node_family_info.name = node_family_info.name;
406 face_connected_node_family_info.item_kind = node_family_info.item_kind;
407 face_connected_node_family_info.connectivity_name =
"FaceToNodes";
408 face_connected_node_family_info.nb_connected_items_per_item = vtk_reader.faceNbNodes();
409 face_connected_node_family_info.connected_items_uids = vtk_reader.faceNodes();
411 auto& face_connected_edge_family_info = face_family_info.connected_family_infos[face_connected_family_index];
412 face_connected_edge_family_info.name = edge_family_info.name;
413 face_connected_edge_family_info.item_kind = edge_family_info.item_kind;
414 face_connected_edge_family_info.connectivity_name =
"FaceToEdges";
415 face_connected_edge_family_info.nb_connected_items_per_item = vtk_reader.faceNbEdges();
416 face_connected_edge_family_info.connected_items_uids = vtk_reader.faceEdges();
418 auto edge_connected_family_index = 0;
419 auto& edge_connected_cell_family_info = edge_family_info.connected_family_infos[edge_connected_family_index++];
420 edge_connected_cell_family_info.name = cell_family_info.name;
421 edge_connected_cell_family_info.item_kind = cell_family_info.item_kind;
422 edge_connected_cell_family_info.connectivity_name =
"EdgeToCells";
423 edge_connected_cell_family_info.nb_connected_items_per_item = vtk_reader.edgeNbCells();
424 edge_connected_cell_family_info.connected_items_uids = vtk_reader.edgeCells();
426 auto& edge_connected_face_family_info = edge_family_info.connected_family_infos[edge_connected_family_index++];
427 edge_connected_face_family_info.name = face_family_info.name;
428 edge_connected_face_family_info.item_kind = face_family_info.item_kind;
429 edge_connected_face_family_info.connectivity_name =
"EdgeToFaces";
430 edge_connected_face_family_info.nb_connected_items_per_item = vtk_reader.edgeNbFaces();
431 edge_connected_face_family_info.connected_items_uids = vtk_reader.edgeFaces();
433 auto& edge_connected_node_family_info = edge_family_info.connected_family_infos[edge_connected_family_index++];
434 edge_connected_node_family_info.name = node_family_info.name;
435 edge_connected_node_family_info.item_kind = node_family_info.item_kind;
436 edge_connected_node_family_info.connectivity_name =
"EdgeToNodes";
437 edge_connected_node_family_info.nb_connected_items_per_item = vtk_reader.edgeNbNodes();
438 edge_connected_node_family_info.connected_items_uids = vtk_reader.edgeNodes();
440 auto node_connected_family_index = 0;
441 auto& node_connected_cell_family_info = node_family_info.connected_family_infos[node_connected_family_index++];
442 node_connected_cell_family_info.name = cell_family_info.name;
443 node_connected_cell_family_info.item_kind = cell_family_info.item_kind;
444 node_connected_cell_family_info.connectivity_name =
"NodeToCells";
445 node_connected_cell_family_info.nb_connected_items_per_item = vtk_reader.nodeNbCells();
446 node_connected_cell_family_info.connected_items_uids = vtk_reader.nodeCells();
448 auto& node_connected_face_family_info = node_family_info.connected_family_infos[node_connected_family_index++];
449 node_connected_face_family_info.name = face_family_info.name;
450 node_connected_face_family_info.item_kind = face_family_info.item_kind;
451 node_connected_face_family_info.connectivity_name =
"NodeToFaces";
452 node_connected_face_family_info.nb_connected_items_per_item = vtk_reader.nodeNbFaces();
453 node_connected_face_family_info.connected_items_uids = vtk_reader.nodeFaces();
455 auto& node_connected_edge_family_info = node_family_info.connected_family_infos[node_connected_family_index++];
456 node_connected_edge_family_info.name = edge_family_info.name;
457 node_connected_edge_family_info.item_kind = edge_family_info.item_kind;
458 node_connected_edge_family_info.connectivity_name =
"NodeToEdges";
459 node_connected_edge_family_info.nb_connected_items_per_item = vtk_reader.nodeNbEdges();
460 node_connected_edge_family_info.connected_items_uids = vtk_reader.nodeEdges();
462 node_family_info.item_coordinates_variable_name =
"NodeCoord";
463 node_family_info.item_coordinates = vtk_reader.nodeCoords();
469void VtkPolyhedralMeshIOService::
473 OStringStream created_infos_str;
474 created_infos_str() <<
"<?xml version='1.0' ?>\n";
475 created_infos_str() <<
"<infos>";
477 if (
auto* cell_data = reader.cellData(); cell_data) {
480 std::iota(vtk_to_arcane_lids.begin(), vtk_to_arcane_lids.end(), 0);
482 for (
auto array_index = 0; array_index < cell_data->GetNumberOfArrays(); ++array_index) {
483 auto* cell_array = cell_data->GetArray(array_index);
486 if (String name = cell_array->GetName(); name.substring(0, 6) ==
"GROUP_") {
487 _createGroup(cell_array, name.substring(6), mesh, mesh->cellFamily(), vtk_to_arcane_lids.constSpan());
488 created_infos_str() <<
"<cell-group name='" << name.substring(6) <<
"'/>";
491 auto var_info = _createVariable(cell_array, name, mesh, mesh->cellFamily(), arcane_to_vtk_lids);
492 created_infos_str() <<
"<cell-variable name='" << name <<
"' "
493 <<
" data-type='" <<
dataTypeName(var_info.m_type) <<
"' "
494 <<
" is_array='" << std::boolalpha << var_info.is_array <<
"'/>";
496 if (m_print_info_level.print_debug_info) {
498 for (
auto tuple_index = 0; tuple_index < cell_array->GetNumberOfTuples(); ++tuple_index) {
499 for (
auto component_index = 0; component_index < cell_array->GetNumberOfComponents(); ++component_index) {
500 debug(
Trace::High) << cell_array->GetName() <<
"[" << tuple_index <<
"][" << component_index <<
"] = " << cell_array->GetComponent(tuple_index, component_index);
507 if (
auto* point_data = reader.pointData(); point_data) {
510 std::iota(vtk_to_arcane_lids.begin(), vtk_to_arcane_lids.end(), 0);
512 for (
auto array_index = 0; array_index < point_data->GetNumberOfArrays(); ++array_index) {
513 auto* point_array = point_data->GetArray(array_index);
514 if (String name = point_array->GetName(); name.substring(0, 6) ==
"GROUP_") {
515 _createGroup(point_array, name.substring(6), mesh, mesh->nodeFamily(), vtk_to_arcane_lids.constSpan());
516 created_infos_str() <<
"<node-group name='" << name.substring(6) <<
"'/>";
519 auto var_info = _createVariable(point_array, name, mesh, mesh->nodeFamily(), arcane_to_vtk_lids);
520 created_infos_str() <<
"<node-variable name='" << name <<
"' "
521 <<
" data-type='" <<
dataTypeName(var_info.m_type) <<
"' "
522 <<
" is_array='" << std::boolalpha << var_info.is_array <<
"'/>";
524 if (m_print_info_level.print_debug_info) {
526 for (
auto tuple_index = 0; tuple_index < point_array->GetNumberOfTuples(); ++tuple_index) {
527 for (
auto component_index = 0; component_index < point_array->GetNumberOfComponents(); ++component_index) {
528 debug(
Trace::High) << point_array->GetName() <<
"[" << tuple_index <<
"][" << component_index <<
"] = " << point_array->GetComponent(tuple_index, component_index);
535 if (
auto* face_data = reader.faceData(); face_data) {
539 _computeFaceVtkArcaneLidConversion(vtk_to_Arcane_lids, arcane_to_vtk_lids, reader, mesh);
540 for (
auto array_index = 0; array_index < face_data->GetNumberOfArrays(); ++array_index) {
541 auto* face_array = face_data->GetArray(array_index);
542 if (String name = face_array->GetName(); name.substring(0, 6) ==
"GROUP_") {
543 _createGroup(face_array, name.substring(6), mesh, mesh->faceFamily(), vtk_to_Arcane_lids);
544 created_infos_str() <<
"<face-group name='" << name.substring(6) <<
"'/>";
547 auto var_info = _createVariable(face_array, name, mesh, mesh->faceFamily(), arcane_to_vtk_lids);
548 created_infos_str() <<
"<face-variable name='" << name <<
"' "
549 <<
" data-type='" <<
dataTypeName(var_info.m_type) <<
"' "
550 <<
" is_array='" << std::boolalpha << var_info.is_array <<
"'/>";
552 if (m_print_info_level.print_debug_info) {
554 for (
auto tuple_index = 0; tuple_index < face_array->GetNumberOfTuples(); ++tuple_index) {
555 for (
auto component_index = 0; component_index < face_array->GetNumberOfComponents(); ++component_index) {
556 debug(
Trace::High) << face_array->GetName() <<
"[" << tuple_index <<
"][" << component_index <<
"] = " << face_array->GetComponent(tuple_index, component_index);
562 created_infos_str() <<
"</infos>";
565 auto* pm = mesh->parallelMng();
566 if (!reader.isEmpty() && pm->isMasterIO()) {
567 String str = created_infos_str.str();
570 bytes.resize(len + 1);
573 pm->broadcastMemoryBuffer(bytes, pm->masterIORank());
574 if (reader.isEmpty()) {
576 XmlNode doc_node = doc->documentNode();
577 _createEmptyVariablesAndGroups(mesh, doc_node);
584void VtkPolyhedralMeshIOService::
590 if (group_items->GetNumberOfComponents() != 1)
591 fatal() << String::format(
"Cannot create item group {0}. Group information in data property must be a scalar", group_name);
592 debug() <<
"Create group " << group_name;
594 arcane_lids.
reserve((
int)group_items->GetNumberOfValues());
595 using GroupDispatcher = vtkArrayDispatch::DispatchByValueType<vtkArrayDispatch::Integrals>;
596 auto group_creator = [&arcane_lids, &vtkToArcaneLid](
auto* array) {
597 vtkIdType numTuples = array->GetNumberOfTuples();
598 vtkDataArrayAccessor<std::remove_pointer_t<
decltype(array)>> array_accessor{ array };
600 for (vtkIdType tupleIdx = 0; tupleIdx < numTuples; ++tupleIdx) {
601 auto value = array_accessor.Get(tupleIdx, 0);
603 arcane_lids.push_back(vtkToArcaneLid[local_id]);
607 if (!GroupDispatcher::Execute(group_items, group_creator))
608 ARCANE_FATAL(
"Cannot create item group {0}. Group information in data property must be an integral type", group_name);
609 debug() <<
" local ids for item group " << group_name <<
" " << arcane_lids;
610 item_family->createGroup(group_name, arcane_lids);
616template <
typename vtkType>
619 using type = vtkType;
634template <
typename T>
using to_arcane_type_t =
typename ToArcaneType<T>::type;
643 if (item_values->GetNumberOfTuples() != item_family->
nbItem())
644 ARCANE_FATAL(
"Cannot create variable {0}, {1} values are given for {2} items in {3} family",
645 variable_name, item_values->GetNumberOfTuples(), item_family->
nbItem(), item_family->
name());
646 debug() <<
"Create mesh variable " << variable_name;
648 auto variable_creator = [
mesh, variable_name, item_family, arcane_to_vtk_lids,
this, &var_info](
auto* values) {
650 using ValueType =
typename std::remove_pointer_t<
decltype(values)>::ValueType;
651 auto* var =
new ItemVariableScalarRefT<to_arcane_type_t<ValueType>>{ vbi, item_family->itemKind() };
653 vtkDataArrayAccessor<std::remove_pointer_t<
decltype(values)>> values_accessor{ values };
655 (*var)[item] = (to_arcane_type_t<ValueType>)values_accessor.Get(arcane_to_vtk_lids[item.localId()], 0);
657 var_info.m_type = DataTypeTraitsT<to_arcane_type_t<ValueType>>::type();
658 var_info.is_array =
false;
660 auto array_variable_creator = [mesh, variable_name, item_family, arcane_to_vtk_lids,
this, &var_info](
auto* values) {
661 VariableBuildInfo vbi{ mesh, variable_name };
662 using ValueType =
typename std::remove_pointer_t<
decltype(values)>::ValueType;
663 auto* var =
new ItemVariableArrayRefT<to_arcane_type_t<ValueType>>{ vbi, item_family->itemKind() };
665 vtkDataArrayAccessor<std::remove_pointer_t<
decltype(values)>> values_accessor{ values };
666 var->resize(values->GetNumberOfComponents());
669 for (
auto& var_value : (*var)[item]) {
670 var_value = (to_arcane_type_t<ValueType>)values_accessor.Get(arcane_to_vtk_lids[item.localId()], index++);
673 var_info.m_type = DataTypeTraitsT<to_arcane_type_t<ValueType>>::type();
674 var_info.is_array =
true;
677 using ValueTypes = vtkTypeList_Create_6(
double,
float,
int,
long,
long long,
short);
678 using ArrayDispatcher = vtkArrayDispatch::DispatchByValueType<ValueTypes>;
680 bool is_variable_created =
false;
681 if (item_values->GetNumberOfComponents() == 1) {
682 is_variable_created = ArrayDispatcher::Execute(item_values, variable_creator);
686 is_variable_created = ArrayDispatcher::Execute(item_values, array_variable_creator);
688 if (!is_variable_created)
689 ARCANE_FATAL(
"Cannot create variable {0}, it's data type is not supported. Only real and integral types are supported", variable_name);
696void VtkPolyhedralMeshIOService::
699 auto face_nodes_unique_ids = reader.faceNodesInFaceMesh();
700 auto face_nb_nodes = reader.faceNbNodesInFaceMesh();
701 auto current_face_index = 0;
702 auto current_face_index_in_face_nodes = 0;
704 mesh->nodeFamily()->itemsUniqueIdToLocalId(face_nodes_local_ids, face_nodes_unique_ids);
705 for (
auto current_face_nb_node : face_nb_nodes) {
706 auto current_face_nodes = face_nodes_local_ids.subConstView(current_face_index_in_face_nodes, current_face_nb_node);
707 current_face_index_in_face_nodes += current_face_nb_node;
708 Node face_first_node{ mesh->nodeFamily()->view()[current_face_nodes[0]] };
709 face_vtk_to_arcane_lids[current_face_index] = MeshUtils::getFaceFromNodesLocalId(face_first_node, current_face_nodes).localId();
710 ++current_face_index;
713 for (
auto arcane_lid : face_vtk_to_arcane_lids) {
714 arcane_to_vtk_lids[arcane_lid] = vtk_lid;
722template <
template <
class>
class VariableRootType>
725 bool has_error =
false;
728 ARCANE_FATAL(
"Invalid data type name {0} for Variable creation in VtkPolyhedralMeshIOService");
730 switch (var_data_type) {
741 ARCANE_FATAL(
"Handle only DT_Int32, DT_Int64, DT_Real in VtkPolyhedralMeshIOService");
748template <
template <
class>
class VariableRootType,
template <
class>
class ArrayVariableRootType>
749void VtkPolyhedralMeshIOService::
752 ARCANE_CHECK_PTR(mesh);
754 for (XmlNode xnode : item_variables_node) {
755 String name = xnode.attrValue(
"name");
756 debug() <<
"Create mesh variable: " << name;
757 String data_type_name = xnode.attrValue(
"data-type");
758 bool is_array = xnode.attrValue(
"is_array") ==
"true";
759 VariableRef* var =
nullptr;
761 var = _createVar<ArrayVariableRootType>(mesh, name, data_type_name, item_kind);
764 var = _createVar<VariableRootType>(mesh, name, data_type_name, item_kind);
766 mesh->variableMng()->_internalApi()->addAutoDestroyVariable(var);
774void VtkPolyhedralMeshIOService::
778 for (XmlNode xnode : groups_node) {
779 String name = xnode.attrValue(
"name");
780 info() <<
"Building group: " << name;
781 item_family->createGroup(name);
788void VtkPolyhedralMeshIOService::
791 ARCANE_CHECK_PTR(mesh);
792 auto document_node = variable_and_group_info.documentElement();
793 _createEmptyVariables<ItemVariableScalarRefT, ItemVariableArrayRefT>(mesh, document_node.children(
"cell-variable"),
IK_Cell);
794 _createEmptyVariables<ItemVariableScalarRefT, ItemVariableArrayRefT>(mesh, document_node.children(
"node-variable"),
IK_Node);
795 _createEmptyVariables<ItemVariableScalarRefT, ItemVariableArrayRefT>(mesh, document_node.children(
"face-variable"),
IK_Face);
797 _createEmptyGroups(mesh, document_node.children(
"cell-group"), mesh->itemFamily(
IK_Cell));
798 _createEmptyGroups(mesh, document_node.children(
"node-group"), mesh->itemFamily(
IK_Node));
799 _createEmptyGroups(mesh, document_node.children(
"face-group"), mesh->itemFamily(
IK_Face));
805VtkPolyhedralMeshIOService::VtkReader::
807: m_filename{ filename }
808, m_print_info_level{ print_info_level }
812 if (filename.empty()) {
813 m_read_status.failure =
true;
814 m_read_status.failure_message =
"filename for polyhedral vtk mesh is empty.";
817 if (filename.endsWith(
"vtk"))
818 _readPlainTextVtkGrid(filename);
819 else if (filename.endsWith(
"vtu"))
820 _readXlmVtkGrid(filename);
822 m_read_status.failure =
true;
823 m_read_status.failure_message = String::format(
"Unsupported vtk extension for file {0}. Supported vtk extension for Polyhedral meshes are {1}",
824 filename, VtkReader::supportedVtkExtensions());
828 m_read_status.failure =
true;
829 m_read_status.failure_message = String::format(
"Cannot read vtk polyhedral file {0}. Vtk grid was not created.", filename);
832 if (m_vtk_grid->GetNumberOfCells() == 0) {
833 m_read_status.failure =
true;
834 m_read_status.failure_message = String::format(
"Cannot read vtk polyhedral file {0}. No cells were found.", filename);
837 if (!m_vtk_grid->GetFaces()) {
838 m_read_status.failure =
true;
839 m_read_status.failure_message = String::format(
"The given mesh vtk file {0} is not a polyhedral mesh, cannot read it", filename);
843 m_cell_data = m_vtk_grid->GetCellData();
844 m_point_data = m_vtk_grid->GetPointData();
847 String faces_filename = m_filename +
"faces.vtk";
848 std::ifstream ifile(faces_filename.localstr());
850 _readPlainTextVtkFaceGrid(faces_filename);
853 faces_filename = m_filename +
"faces.vtp";
854 ifile = std::ifstream{ faces_filename.localstr() };
856 _readXmlVtkFaceGrid(faces_filename);
861 faces_filename.split(faces_filename_and_extension,
'.');
864 m_read_status.info_message = String::format(
"Information no face mesh given {0}{1} (.vtk or .vtp) to define face variables or groups on faces.",
865 faces_filename_and_extension[0],
866 faces_filename_and_extension[1]);
869 if (m_vtk_face_grid) {
870 if (m_vtk_face_grid->GetNumberOfCells() == 0) {
871 m_read_status.failure =
true;
872 m_read_status.failure_message = m_read_status.failure_message + String::format(
" Error in reading face information for groups in mesh file {0} ", faces_filename);
875 m_face_data = m_vtk_face_grid->GetCellData();
876 m_poly_data = m_vtk_face_grid->GetPolys();
880 m_read_status.failure =
true;
881 m_read_status.failure_message = m_read_status.failure_message + String::format(
"Face data could not be built from file {0}{1} (.vtk or .vtu).", faces_filename_and_extension[0], faces_filename_and_extension[1]);
885 if (m_print_info_level.print_mesh_info)
897 if (m_cell_uids.empty()) {
899 m_cell_uids.reserve(m_vtk_grid->GetNumberOfCells());
900 m_cell_nb_nodes.reserve(m_vtk_grid->GetNumberOfCells());
901 m_cell_node_uids.reserve(10 * m_vtk_grid->GetNumberOfCells());
902 auto* cell_iter = m_vtk_grid->NewCellIterator();
903 cell_iter->InitTraversal();
904 while (!cell_iter->IsDoneWithTraversal()) {
905 m_cell_uids.push_back(cell_iter->GetCellId());
906 m_cell_nb_nodes.push_back(
Integer(cell_iter->GetNumberOfPoints()));
907 ArrayView<vtkIdType> cell_nodes{
Integer(cell_iter->GetNumberOfPoints()), cell_iter->GetPointIds()->GetPointer(0) };
908 std::for_each(cell_nodes.begin(), cell_nodes.end(), [
this](
auto uid) { this->m_cell_node_uids.push_back(uid); });
909 cell_iter->GoToNextCell();
923 if (m_node_uids.empty()) {
925 auto nb_nodes = m_vtk_grid->GetNumberOfPoints();
926 m_node_uids.resize(nb_nodes);
927 m_node_nb_cells.resize(nb_nodes);
928 m_node_cell_uids.reserve(8 * nb_nodes);
929 m_node_uid_to_index.resize(nb_nodes);
930 for (
int node_index = 0; node_index < nb_nodes; ++node_index) {
931 Int64 node_uid = node_index;
932 m_node_uids[node_index] = node_uid;
933 auto cell_nodes = vtkIdList::New();
934 m_vtk_grid->GetPointCells(node_index, cell_nodes);
935 Int64Span cell_nodes_view((
Int64*)cell_nodes->GetPointer(0), cell_nodes->GetNumberOfIds());
936 m_node_cell_uids.addRange(cell_nodes_view);
937 m_node_nb_cells[node_index] = (
Int32)cell_nodes->GetNumberOfIds();
940 m_node_uid_to_index[node_uid] = node_index;
957 if (!m_face_uids.empty())
960 auto* cell_iter = m_vtk_grid->NewCellIterator();
961 cell_iter->InitTraversal();
962 vtkIdType nb_face_estimation = 0;
963 while (!cell_iter->IsDoneWithTraversal()) {
964 vtkIdType cell_nb_faces = 0;
965 vtkIdType_generic* points{
nullptr };
966 m_vtk_grid->GetFaceStream(cell_iter->GetCellId(), cell_nb_faces, points);
967 nb_face_estimation += cell_nb_faces;
968 cell_iter->GoToNextCell();
970 m_face_uids.reserve(nb_face_estimation);
971 auto const* faces = m_vtk_grid->GetFaces();
976 ARCANE_FATAL(
"Mesh {0} is not polyhedral: faces are not defined", m_filename);
979 auto face_info_size = faces->GetNumberOfValues();
980 m_face_node_uids.reserve(face_info_size);
981 m_face_nb_nodes.reserve(nb_face_estimation);
982 m_face_cell_uids.reserve(2 * nb_face_estimation);
983 m_face_nb_cells.reserve(nb_face_estimation);
984 m_cell_face_uids.reserve(8 * m_cell_uids.size());
985 m_cell_nb_faces.resize(m_cell_uids.size(), 0);
986 m_cell_face_indexes.resize(m_cell_uids.size(), -1);
987 m_face_uid_indexes.resize(2 * nb_face_estimation, -1);
989 current_face_nodes.
reserve(10);
990 sorted_current_face_nodes.reserve(10);
991 UniqueArray<UniqueArray<Int64>> node_faces(m_node_uids.size());
992 UniqueArray<Int32> face_offsets;
993 face_offsets.reserve(nb_face_estimation);
994 face_offsets.push_back(0);
995 FaceUidToIndexMap face_uid_to_index;
996 face_uid_to_index.reserve(nb_face_estimation);
998 auto cell_face_index = 0;
999 auto global_face_index = 0;
1000 auto face_uid_index = 0;
1001 for (
int face_info_index = 0; face_info_index < face_info_size; cell_index++) {
1002 auto current_cell_nb_faces =
Int32(faces->GetValue(face_info_index++));
1003 m_cell_face_indexes[m_cell_uids[cell_index]] = cell_face_index;
1004 for (
auto face_index = 0; face_index < current_cell_nb_faces; ++face_index, ++global_face_index) {
1005 auto current_face_nb_nodes =
Int32(faces->GetValue(face_info_index++));
1006 m_cell_nb_faces[m_cell_uids[cell_index]] += 1;
1007 for (
int node_index = 0; node_index < current_face_nb_nodes; ++node_index) {
1008 current_face_nodes.push_back(faces->GetValue(face_info_index++));
1010 sorted_current_face_nodes.resize(current_face_nodes.size());
1011 auto is_front_cell = mesh_utils::reorderNodesOfFace(current_face_nodes, sorted_current_face_nodes);
1012 auto [is_face_found, existing_face_index] = _findFace(sorted_current_face_nodes, node_faces,
1013 m_node_uid_to_index, m_face_nb_nodes,
1014 face_uid_to_index, face_offsets, m_face_node_uids);
1015 if (!is_face_found) {
1016 for (
auto node_uid : current_face_nodes) {
1017 node_faces[m_node_uid_to_index[node_uid]].push_back(face_uid);
1019 m_cell_face_uids.push_back(face_uid);
1020 m_face_uids.push_back(face_uid);
1021 m_face_nb_nodes.push_back(current_face_nb_nodes);
1022 m_face_node_uids.addRange(sorted_current_face_nodes);
1023 m_face_nb_cells.push_back(1);
1024 m_face_uid_indexes[global_face_index] = face_uid_index;
1025 face_uid_to_index.push_back(face_uid_index);
1026 auto previous_offset = face_offsets.back();
1027 face_offsets.push_back(previous_offset + sorted_current_face_nodes.size());
1030 if (is_front_cell) {
1031 m_face_cell_uids.push_back(NULL_ITEM_UNIQUE_ID);
1032 m_face_cell_uids.push_back(m_cell_uids[cell_index]);
1035 m_face_cell_uids.push_back(m_cell_uids[cell_index]);
1036 m_face_cell_uids.push_back(NULL_ITEM_UNIQUE_ID);
1040 m_cell_face_uids.push_back(m_face_uids[existing_face_index]);
1041 m_face_nb_cells[existing_face_index] += 1;
1042 m_face_uid_indexes[global_face_index] = existing_face_index;
1044 if (is_front_cell) {
1045 if (m_face_cell_uids[2 * existing_face_index + 1] != NULL_ITEM_UNIQUE_ID) {
1046 ARCANE_FATAL(
"Problem in face orientation, face uid {0}, nodes {1}, same orientation in cell {2} and {3}. Change mesh file.",
1047 m_face_uids[existing_face_index],
1049 m_face_cell_uids[2 * existing_face_index + 1],
1050 m_cell_uids[cell_index]);
1052 m_face_cell_uids[2 * existing_face_index + 1] = m_cell_uids[cell_index];
1055 if (m_face_cell_uids[2 * existing_face_index] != NULL_ITEM_UNIQUE_ID) {
1056 ARCANE_FATAL(
"Problem in face orientation, face uid {0}, nodes {1}, same orientation in cell {2} and {3}. Change mesh file.",
1057 m_face_uids[existing_face_index],
1059 m_face_cell_uids[2 * existing_face_index],
1060 m_cell_uids[cell_index]);
1062 m_face_cell_uids[2 * existing_face_index] = m_cell_uids[cell_index];
1065 current_face_nodes.clear();
1066 sorted_current_face_nodes.clear();
1068 cell_face_index += m_cell_nb_faces[m_cell_uids[cell_index]];
1071 m_node_nb_faces.resize(m_node_uids.size(), 0);
1072 _flattenConnectivity(node_faces.constSpan(), m_node_nb_faces, m_node_face_uids);
1074 if (m_print_info_level.print_debug_info) {
1075 std::cout <<
"================FACE NODES ==============" << std::endl;
1076 std::copy(m_face_node_uids.begin(), m_face_node_uids.end(), std::ostream_iterator<Int64>(std::cout,
" "));
1077 std::cout << std::endl;
1078 std::copy(m_face_nb_nodes.begin(), m_face_nb_nodes.end(), std::ostream_iterator<Int64>(std::cout,
" "));
1079 std::cout << std::endl;
1080 std::copy(m_cell_face_indexes.begin(), m_cell_face_indexes.end(), std::ostream_iterator<Int64>(std::cout,
" "));
1081 std::cout << std::endl;
1095 if (!m_edge_uids.empty())
1101 m_edge_uids.reserve(2 * m_vtk_grid->GetNumberOfPoints());
1102 auto const* faces = m_vtk_grid->GetFaces();
1107 ARCANE_FATAL(
"Mesh {0} is not polyhedral: faces are not defined", m_filename);
1110 auto nb_edge_estimation = 2 * m_edge_uids.capacity();
1111 m_edge_node_uids.reserve(nb_edge_estimation);
1112 auto face_info_size = faces->GetNumberOfValues();
1113 auto cell_index = 0;
1114 auto global_face_index = 0;
1115 auto new_edge_index = 0;
1116 UniqueArray<std::set<Int64>> edge_cells;
1117 UniqueArray<Int64UniqueArray> edge_faces;
1118 edge_cells.reserve(m_edge_uids.capacity());
1119 edge_faces.reserve(m_edge_uids.capacity());
1120 m_cell_nb_edges.resize(m_cell_uids.size(), 0);
1121 m_cell_edge_uids.reserve(20 * m_cell_uids.size());
1122 UniqueArray<std::set<Int64>> face_edges;
1123 face_edges.resize(m_face_uids.size());
1124 UniqueArray<std::set<Int64>> cell_edges;
1125 cell_edges.resize(m_cell_uids.size());
1126 UniqueArray<Int64UniqueArray> node_edges;
1127 node_edges.resize(m_node_uids.size());
1128 EdgeUidToIndexMap edge_uid_to_index;
1129 edge_uid_to_index.reserve(nb_edge_estimation);
1131 edge_offsets.
reserve(nb_edge_estimation);
1132 edge_offsets.push_back(0);
1133 m_edge_nb_nodes.reserve(nb_edge_estimation);
1134 for (
int face_info_index = 0; face_info_index < face_info_size; ++cell_index) {
1135 auto current_cell_nb_faces =
Int32(faces->GetValue(face_info_index++));
1136 for (
auto face_index = 0; face_index < current_cell_nb_faces; ++face_index, ++global_face_index) {
1137 auto current_face_nb_nodes =
Int32(faces->GetValue(face_info_index++));
1138 auto first_face_node_uid =
Int32(faces->GetValue(face_info_index));
1139 UniqueArray<Int64> current_edge(2), sorted_edge(2);
1140 for (
int node_index = 0; node_index < current_face_nb_nodes - 1; ++node_index) {
1141 current_edge = UniqueArray<Int64>{ faces->GetValue(face_info_index++), faces->GetValue(face_info_index) };
1142 mesh_utils::reorderNodesOfFace(current_edge, sorted_edge);
1143 auto [is_edge_found, existing_edge_index] = _findFace(sorted_edge, node_edges,
1144 m_node_uid_to_index,
1146 edge_uid_to_index, edge_offsets, m_edge_node_uids);
1147 if (!is_edge_found) {
1148 m_cell_nb_edges[cell_index] += 1;
1149 face_edges[m_face_uid_indexes[global_face_index]].insert(edge_uid);
1150 cell_edges[cell_index].insert(edge_uid);
1151 for (
auto node : current_edge) {
1152 node_edges[node].push_back(edge_uid);
1154 edge_cells.push_back(std::set{ m_cell_uids[cell_index] });
1155 edge_faces.push_back(
Int64UniqueArray{ m_cell_face_uids[m_cell_face_indexes[cell_index] + face_index] });
1156 m_edge_uids.push_back(edge_uid++);
1157 m_edge_node_uids.addRange(sorted_edge);
1158 edge_uid_to_index.push_back(new_edge_index);
1159 auto current_offset = edge_offsets.back();
1160 edge_offsets.push_back(current_offset + 2);
1161 m_edge_nb_nodes.push_back(2);
1165 edge_cells[existing_edge_index].insert(m_cell_uids[cell_index]);
1166 edge_faces[existing_edge_index].push_back(m_cell_face_uids[m_cell_face_indexes[cell_index] + face_index]);
1167 face_edges[m_face_uid_indexes[global_face_index]].insert(m_edge_uids[existing_edge_index]);
1168 cell_edges[cell_index].insert(m_edge_uids[existing_edge_index]);
1171 current_edge = UniqueArray<Int64>{ faces->GetValue(face_info_index++), first_face_node_uid };
1172 mesh_utils::reorderNodesOfFace(current_edge, sorted_edge);
1173 auto [is_edge_found, existing_edge_index] = _findFace(sorted_edge, node_edges,
1174 m_node_uid_to_index,
1176 edge_uid_to_index, edge_offsets, m_edge_node_uids);
1177 if (!is_edge_found) {
1178 m_cell_nb_edges[cell_index] += 1;
1179 edge_cells.push_back(std::set{ m_cell_uids[cell_index] });
1180 edge_faces.push_back(
Int64UniqueArray{ m_cell_face_uids[m_cell_face_indexes[cell_index] + face_index] });
1181 face_edges[m_face_uid_indexes[global_face_index]].insert(edge_uid);
1182 cell_edges[cell_index].insert(edge_uid);
1183 for (
auto node : current_edge) {
1184 node_edges[node].push_back(edge_uid);
1186 m_edge_uids.push_back(edge_uid++);
1187 m_edge_node_uids.addRange(sorted_edge);
1188 edge_uid_to_index.push_back(new_edge_index);
1189 auto current_offset = edge_offsets.back();
1190 edge_offsets.push_back(current_offset + 2);
1191 m_edge_nb_nodes.push_back(2);
1195 edge_cells[existing_edge_index].insert(m_cell_uids[cell_index]);
1196 edge_faces[existing_edge_index].push_back(m_cell_face_uids[m_cell_face_indexes[cell_index] + face_index]);
1197 face_edges[m_face_uid_indexes[global_face_index]].insert(m_edge_uids[existing_edge_index]);
1198 cell_edges[cell_index].insert(m_edge_uids[existing_edge_index]);
1203 m_edge_nb_cells.resize(m_edge_uids.size(), 0);
1204 _flattenConnectivity(edge_cells.constSpan(), m_edge_nb_cells, m_edge_cell_uids);
1207 m_edge_nb_faces.resize(m_edge_uids.size(), 0);
1208 _flattenConnectivity(edge_faces.constSpan(), m_edge_nb_faces, m_edge_face_uids);
1211 m_face_nb_edges.resize(m_face_uids.size(), 0);
1212 _flattenConnectivity(face_edges.constSpan(), m_face_nb_edges, m_face_edge_uids);
1215 m_cell_nb_edges.resize(m_cell_uids.size(), 0);
1216 _flattenConnectivity(cell_edges, m_cell_nb_edges, m_cell_edge_uids);
1219 m_node_nb_edges.resize(m_node_uids.size(), 0);
1220 _flattenConnectivity(node_edges, m_node_nb_edges, m_node_edge_uids);
1222 if (m_print_info_level.print_debug_info) {
1223 std::cout <<
"================EDGE NODES ==============" << std::endl;
1224 std::copy(m_edge_node_uids.begin(), m_edge_node_uids.end(), std::ostream_iterator<Int64>(std::cout,
" "));
1225 std::cout << std::endl;
1226 std::cout <<
"================FACE EDGES ==============" << std::endl;
1227 std::copy(m_face_nb_edges.begin(), m_face_nb_edges.end(), std::ostream_iterator<Int32>(std::cout,
" "));
1228 std::cout << std::endl;
1229 std::copy(m_face_edge_uids.begin(), m_face_edge_uids.end(), std::ostream_iterator<Int64>(std::cout,
" "));
1230 std::cout << std::endl;
1231 std::cout <<
"================CELL EDGES ==============" << std::endl;
1232 std::copy(m_cell_nb_edges.begin(), m_cell_nb_edges.end(), std::ostream_iterator<Int32>(std::cout,
" "));
1233 std::cout << std::endl;
1234 std::copy(m_cell_edge_uids.begin(), m_cell_edge_uids.end(), std::ostream_iterator<Int64>(std::cout,
" "));
1235 std::cout << std::endl;
1243std::pair<bool, Int32> VtkPolyhedralMeshIOService::VtkReader::
1245 const UniqueArray<Int64UniqueArray>& node_face_uids,
1246 const NodeUidToIndexMap& node_uid_to_index,
1248 const FaceUidToIndexMap& face_uid_to_index,
1249 const UniqueArray<Int32>& face_offsets,
1252 auto first_node_uid = sorted_face_nodes[0];
1253 auto first_node_index = node_uid_to_index[first_node_uid];
1255 for (
auto face_uid : node_face_uids[first_node_index]) {
1256 auto face_index = face_uid_to_index[face_uid];
1257 auto face_offset = face_offsets[face_index];
1258 auto face_nb_node = face_nb_nodes[face_index];
1259 if (face_nb_node == sorted_face_nodes.size()) {
1260 bool is_same_face =
true;
1261 for (
auto index = 0; index < face_nb_node; ++index) {
1262 if (sorted_face_nodes[index] != face_node_uids[face_offset + index]) {
1263 is_same_face =
false;
1267 return {
true, face_index };
1270 return {
false, -1 };
1276Integer VtkPolyhedralMeshIOService::VtkReader::
1279 if (m_node_uids.empty())
1281 return m_node_uids.size();
1290 if (m_cell_node_uids.empty())
1292 return m_cell_node_uids;
1301 if (m_cell_nb_nodes.empty())
1303 return m_cell_nb_nodes;
1312 if (m_face_node_uids.empty())
1314 return m_face_node_uids;
1323 if (m_face_nb_nodes.empty())
1325 return m_face_nb_nodes;
1332faceNodesInFaceMesh()
1334 if (m_face_node_uids_in_face_mesh.empty())
1335 _readfaceNodesInFaceMesh();
1336 return m_face_node_uids_in_face_mesh;
1343faceNbNodesInFaceMesh()
1345 if (m_face_nb_nodes_in_face_mesh.empty())
1346 _readfaceNodesInFaceMesh();
1347 return m_face_nb_nodes_in_face_mesh;
1356 if (m_edge_node_uids.empty())
1358 return m_edge_nb_nodes;
1367 if (m_edge_node_uids.empty())
1369 return m_edge_node_uids;
1378 if (m_face_cell_uids.empty())
1381 if (m_print_info_level.print_debug_info) {
1382 std::cout <<
"=================FACE CELLS================="
1384 std::copy(m_face_cell_uids.begin(), m_face_cell_uids.end(), std::ostream_iterator<Int64>(std::cout,
" "));
1386 std::cout <<
"=================END FACE CELLS================="
1389 return m_face_cell_uids;
1398 if (m_face_nb_cells.empty())
1400 return m_face_nb_cells;
1409 if (m_edge_nb_cells.empty())
1411 return m_edge_nb_cells;
1420 if (m_edge_cell_uids.empty())
1422 return m_edge_cell_uids;
1431 if (m_cell_nb_faces.empty())
1433 return m_cell_nb_faces;
1442 if (m_cell_face_uids.empty())
1444 return m_cell_face_uids;
1453 if (m_edge_nb_faces.empty())
1455 return m_edge_nb_faces;
1464 if (m_edge_face_uids.empty())
1466 return m_edge_face_uids;
1475 if (m_cell_nb_edges.empty())
1477 return m_cell_nb_edges;
1486 if (m_cell_edge_uids.empty())
1488 return m_cell_edge_uids;
1497 if (m_face_nb_edges.empty())
1499 return m_face_nb_edges;
1508 if (m_face_edge_uids.empty())
1510 return m_face_edge_uids;
1516template <
typename Connectivity2DArray>
1517void VtkPolyhedralMeshIOService::VtkReader::
1518_flattenConnectivity(Connectivity2DArray connected_item_2darray,
1519 Int32Span nb_connected_item_per_source_item,
1523 std::transform(connected_item_2darray.begin(), connected_item_2darray.end(), nb_connected_item_per_source_item.begin(), [](
auto const& connected_items) {
1524 return connected_items.size();
1527 connected_item_array.reserve(std::accumulate(nb_connected_item_per_source_item.begin(), nb_connected_item_per_source_item.end(), 0));
1528 std::for_each(connected_item_2darray.begin(), connected_item_2darray.end(), [&connected_item_array](
auto const& connected_items) {
1529 for (auto const& connected_item : connected_items) {
1530 connected_item_array.push_back(connected_item);
1541 if (m_node_nb_cells.empty())
1543 return m_node_nb_cells;
1552 if (m_node_cell_uids.empty())
1554 return m_node_cell_uids;
1563 if (m_node_nb_faces.empty())
1565 return m_node_nb_faces;
1574 if (m_node_face_uids.empty())
1576 return m_node_face_uids;
1585 if (m_node_nb_edges.empty())
1587 return m_node_nb_edges;
1596 if (m_node_edge_uids.empty())
1598 return m_node_edge_uids;
1608 return m_node_coordinates;
1609 if (m_node_coordinates.empty()) {
1611 auto point_coords = m_vtk_grid->GetPoints()->GetData();
1612 if (m_print_info_level.print_debug_info) {
1613 std::cout <<
"======= Point COORDS ====" << std::endl;
1614 std::ostringstream oss;
1615 point_coords->PrintSelf(oss, vtkIndent{ 2 });
1616 std::cout << oss.str() << std::endl;
1618 auto nb_nodes = m_vtk_grid->GetNumberOfPoints();
1619 for (
int i = 0; i < nb_nodes; ++i) {
1620 if (m_print_info_level.print_debug_info) {
1621 std::cout <<
"==========current point coordinates : ( ";
1622 std::cout << *(point_coords->GetTuple(i)) <<
" , ";
1623 std::cout << *(point_coords->GetTuple(i) + 1) <<
" , ";
1624 std::cout << *(point_coords->GetTuple(i) + 2) <<
" ) ===" << std::endl;
1626 m_node_coordinates.add({ *(point_coords->GetTuple(i)),
1627 *(point_coords->GetTuple(i) + 1),
1628 *(point_coords->GetTuple(i) + 2) });
1631 return m_node_coordinates;
1637vtkCellData* VtkPolyhedralMeshIOService::VtkReader::
1646vtkPointData* VtkPolyhedralMeshIOService::VtkReader::
1649 return m_point_data;
1655void VtkPolyhedralMeshIOService::VtkReader::
1656_printMeshInfos()
const
1659 std::cout <<
"-- VTK GRID READ "
1660 <<
" NB CELLS " << m_vtk_grid->GetNumberOfCells() << std::endl;
1662 auto* cell_iter = m_vtk_grid->vtkDataSet::NewCellIterator();
1663 cell_iter->InitTraversal();
1664 vtkIdType_generic* cell_faces{
nullptr };
1665 vtkIdType nb_faces = 0;
1666 while (!cell_iter->IsDoneWithTraversal()) {
1667 std::cout <<
"---- visiting cell id " << cell_iter->GetCellId() << std::endl;
1668 std::cout <<
"---- cell number of faces " << cell_iter->GetNumberOfFaces() << std::endl;
1669 std::cout <<
"---- cell number of points " << cell_iter->GetNumberOfPoints() << std::endl;
1670 m_vtk_grid->GetFaceStream(cell_iter->GetCellId(), nb_faces, cell_faces);
1671 for (
auto iface = 0; iface < nb_faces; ++iface) {
1672 auto face_nb_nodes = *cell_faces++;
1673 std::cout <<
"---- has face with " << face_nb_nodes <<
" nodes. Node ids : ";
1674 for (
int inode = 0; inode < face_nb_nodes; ++inode) {
1675 std::cout << *cell_faces++ <<
" ";
1677 std::cout << std::endl;
1679 cell_iter->GoToNextCell();
1686vtkCellData* VtkPolyhedralMeshIOService::VtkReader::faceData()
1696void VtkPolyhedralMeshIOService::VtkReader::
1697_readPlainTextVtkGrid(
const String& filename)
1699 m_vtk_grid_reader->SetFileName(filename.localstr());
1700 m_vtk_grid_reader->ReadAllScalarsOn();
1701 m_vtk_grid_reader->Update();
1702 m_vtk_grid = m_vtk_grid_reader->GetOutput();
1708void VtkPolyhedralMeshIOService::VtkReader::
1709_readXlmVtkGrid(
const String& filename)
1711 m_vtk_xml_grid_reader->SetFileName(filename.localstr());
1712 m_vtk_xml_grid_reader->Update();
1713 m_vtk_grid = m_vtk_xml_grid_reader->GetOutput();
1719void VtkPolyhedralMeshIOService::VtkReader::
1720_readPlainTextVtkFaceGrid(
const String& faces_filename)
1722 m_vtk_face_grid_reader->SetFileName(faces_filename.localstr());
1723 m_vtk_face_grid_reader->ReadAllScalarsOn();
1724 m_vtk_face_grid_reader->Update();
1725 m_vtk_face_grid = m_vtk_face_grid_reader->GetOutput();
1731void VtkPolyhedralMeshIOService::VtkReader::
1732_readXmlVtkFaceGrid(
const String& faces_filename)
1734 m_vtk_xml_face_grid_reader->SetFileName(faces_filename.localstr());
1735 m_vtk_xml_face_grid_reader->Update();
1736 m_vtk_face_grid = m_vtk_xml_face_grid_reader->GetOutput();
1742void VtkPolyhedralMeshIOService::VtkReader::
1743_readfaceNodesInFaceMesh()
1745 m_face_nb_nodes_in_face_mesh.resize(m_poly_data->GetNumberOfCells());
1746 m_face_node_uids_in_face_mesh.reserve(m_poly_data->GetNumberOfCells() * m_poly_data->GetMaxCellSize());
1747 m_poly_data->InitTraversal();
1748 vtkIdType face_nb_nodes;
1749 vtkIdType_generic* face_nodes;
1751 auto face_nb_node_index = 0;
1754 current_face_node_uids.
reserve(m_poly_data->GetMaxCellSize());
1755 reordered_current_face_node_uids.reserve(m_poly_data->GetMaxCellSize());
1757 while (m_poly_data->GetNextCell(face_nb_nodes, face_nodes)) {
1758 m_face_nb_nodes_in_face_mesh[face_nb_node_index] = face_nb_nodes;
1759 ConstArrayView<vtkIdType> face_nodes_view(face_nb_nodes, face_nodes);
1760 current_face_node_uids.resize(face_nb_nodes);
1761 reordered_current_face_node_uids.resize(face_nb_nodes);
1762 std::copy(face_nodes_view.begin(), face_nodes_view.end(), current_face_node_uids.begin());
1763 MeshUtils::reorderNodesOfFace(current_face_node_uids, reordered_current_face_node_uids);
1764 std::copy(reordered_current_face_node_uids.begin(), reordered_current_face_node_uids.end(), std::back_inserter(m_face_node_uids_in_face_mesh));
1765 ++face_nb_node_index;
1771void VtkPolyhedralMeshIOService::VtkReader::
1772_checkVtkGrid()
const
1775 ARCANE_FATAL(
"Polyhedral vtk grid not loaded. Cannot continue.");
Arcane configuration file.
#define ARCANE_CHECK_POINTER(ptr)
Macro returning the pointer ptr if it is not null or throwing an exception if it is null.
#define ARCANE_FATAL(...)
Macro throwing a FatalErrorException.
Declarations of Arcane's general types.
Utility functions for the mesh.
This file contains the various service factories and macros for registering services.
#define ARCANE_SERVICE_INTERFACE(ainterface)
Macro to declare an interface when registering a service.
Management of references to a C++ class.
Generation de la classe de base du Service.
CaseOptionsVtkPolyhedralMeshIO * options() const
Options du jeu de données du service.
ArcaneVtkPolyhedralMeshIOObject(const Arcane::ServiceBuildInfo &sbi)
Constructeur.
void reserve(Int64 new_capacity)
Reserves memory for new_capacity elements.
Necessary information for reading a mesh file.
Interface of an entity family.
virtual String name() const =0
Family name.
virtual Integer nbItem() const =0
Number of entities.
virtual String name() const =0
Mesh name.
Interface of a mesh creation/reading service.
virtual IVariableMng * variableMng() const =0
Associated variable manager.
virtual void addAutoDestroyVariable(VariableRef *var)=0
Adds the variable to the list of variables that are kept until the end of execution.
virtual IVariableMngInternal * _internalApi()=0
Internal Arcane API.
static IXmlDocumentHolder * loadFromBuffer(Span< const Byte > buffer, const String &name, ITraceMng *tm)
Loads an XML document.
Parameters necessary for building a mesh.
MeshBuildInfo & addNeedPartitioning(bool v)
Indicates whether the generator needs to call a partitioner.
MeshBuildInfo & addMeshKind(const MeshKind &v)
Sets the mesh characteristics.
MeshBuildInfo & addFactoryName(const String &factory_name)
Sets the factory name to create this mesh.
Characteristics of a mesh.
Reference to an instance.
Structure containing the information to create a service.
Unicode character string.
TraceAccessor(ITraceMng *m)
Constructs an accessor via the trace manager m.
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flow for a debug message.
TraceMessage fatal() const
Flow for a fatal error message.
TraceMessage info() const
Flow for an information message.
ITraceMng * traceMng() const
Trace manager.
1D data vector with value semantics (STL style).
Parameters necessary for building a variable.
Information characterizing a variable.
void fillMeshBuildInfo(MeshBuildInfo &build_info) override
Fills build_info with the necessary information to create the mesh.
void allocateMeshItems(IPrimaryMesh *pm) override
Allocates the mesh entities managed by this service.
Ref< IMeshBuilder > createBuilder(const CaseMeshReaderReadInfo &read_info) const override
Returns a builder to create and read the mesh whose information is specified in read_info.
List of nodes of a DOM tree.
const value_type & const_reference
Type constant reference of an element in the array.
#define ARCANE_REGISTER_SERVICE(aclass, a_service_property,...)
Macro for registering a service.
Integer len(const char *s)
Returns the length of the string s.
-- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature --
Span< Int64 > Int64Span
C equivalent of a 1D array of 64-bit integers.
UniqueArray< Int64 > Int64UniqueArray
Dynamic 1D array of 64-bit integers.
std::int64_t Int64
Signed integer type of 64 bits.
ArrayView< Real3 > Real3ArrayView
C equivalent of a 1D array of Real3.
Int32 Integer
Type representing an integer.
UniqueArray< Real3 > Real3UniqueArray
Dynamic 1D array of rank 3 vectors.
ConstArrayView< Int32 > Int32ConstArrayView
C equivalent of a 1D array of 32-bit integers.
@ Polyhedral
Polyhedral mesh.
@ ST_CaseOption
The service is used at the dataset level.
ConstArrayView< Int64 > Int64ConstArrayView
C equivalent of a 1D array of 64-bit integers.
UniqueArray< Byte > ByteUniqueArray
Dynamic 1D array of characters.
SharedArray< Int32 > Int32SharedArray
Dynamic 1D array of 32-bit integers.
UniqueArray< Int32 > Int32UniqueArray
Dynamic 1D array of 32-bit integers.
eItemKind
Mesh entity type.
@ IK_Node
Node mesh entity.
@ IK_Cell
Cell mesh entity.
@ IK_Face
Face mesh entity.
@ IK_Edge
Edge mesh entity.
double Real
Type representing a real number.
ConstArrayView< Byte > ByteConstArrayView
C equivalent of a 1D array of characters.
auto makeRef(InstanceType *t) -> Ref< InstanceType >
Creates a reference on a pointer.
UniqueArray< String > StringUniqueArray
Dynamic 1D array of strings.
Span< Int32 > Int32Span
C equivalent of a 1D array of 32-bit integers.
@ DT_Int32
32-bit integer data type
@ DT_Int64
64-bit integer data type
@ DT_Unknown
Unknown or uninitialized data type.
Span< const Int32 > Int32ConstSpan
Read-only view of a 1D array of 32-bit integers.
ARCANE_DATATYPE_EXPORT eDataType dataTypeFromName(const char *name, bool &has_error)
Finds the type associated with name.
const char * dataTypeName(eDataType type)
Data type name.
std::int32_t Int32
Signed integer type of 32 bits.