14#include "arcane/utils/ArcanePrecomp.h"
15#include "arcane/utils/Iostream.h"
16#include "arcane/utils/StdHeader.h"
17#include "arcane/utils/HashTableMap.h"
18#include "arcane/utils/ValueConvert.h"
19#include "arcane/utils/ScopedPtr.h"
20#include "arcane/utils/ITraceMng.h"
21#include "arcane/utils/String.h"
22#include "arcane/utils/IOException.h"
23#include "arcane/utils/Collection.h"
24#include "arcane/utils/Enumerator.h"
25#include "arcane/utils/NotImplementedException.h"
26#include "arcane/utils/Real3.h"
28#include "arcane/core/AbstractService.h"
29#include "arcane/core/FactoryService.h"
30#include "arcane/core/IMainFactory.h"
31#include "arcane/core/IMeshReader.h"
32#include "arcane/core/ISubDomain.h"
33#include "arcane/core/IMesh.h"
34#include "arcane/core/IMeshSubMeshTransition.h"
35#include "arcane/core/IItemFamily.h"
36#include "arcane/core/Item.h"
38#include "arcane/core/VariableTypes.h"
39#include "arcane/core/IVariableAccessor.h"
40#include "arcane/core/IParallelMng.h"
41#include "arcane/core/IIOMng.h"
42#include "arcane/core/IXmlDocumentHolder.h"
43#include "arcane/core/XmlNodeList.h"
44#include "arcane/core/XmlNode.h"
45#include "arcane/core/IMeshUtilities.h"
46#include "arcane/core/IMeshWriter.h"
47#include "arcane/core/BasicService.h"
48#include "arcane/core/ICaseMeshReader.h"
49#include "arcane/core/IMeshBuilder.h"
50#include "arcane/core/ItemPrinter.h"
51#include "arcane/core/MeshKind.h"
52#include "arcane/core/ServiceBuilder.h"
54#include <vtkXMLUnstructuredGridReader.h>
55#include <vtkUnstructuredGrid.h>
57#include <vtkDataObject.h>
58#include <vtkDataSet.h>
59#include <vtkCellData.h>
60#include <vtkPointData.h>
61#include <vtkLongArray.h>
62#include <vtkIntArray.h>
67#define VTK_FILE_EXT_VTI "vti"
68#define VTK_FILE_EXT_VTP "vtp"
69#define VTK_FILE_EXT_VTR "vtr"
70#define VTK_FILE_EXT_VTS "vts"
71#define VTK_FILE_EXT_VTU "vtu"
72#define VTK_FILE_EXT_PVTI "pvti"
73#define VTK_FILE_EXT_PVTP "pvtp"
74#define VTK_FILE_EXT_PVTR "pvtr"
75#define VTK_FILE_EXT_PVTS "pvts"
76#define VTK_FILE_EXT_PVTU "pvtu"
78#if VTK_MAJOR_VERSION >= 9
79#define CURRENT_VTK_VERSION_LONG_TYPE VTK_LONG_LONG
80#include <vtkLongLongArray.h>
81using vtkLongArrayType = vtkLongLongArray;
83#define CURRENT_VTK_VERSION_LONG_TYPE VTK_LONG
84using vtkLongArrayType = vtkLongArray;
99class VtuMeshReaderBase
103 explicit VtuMeshReaderBase(
ITraceMng* trace_mng);
104 virtual ~VtuMeshReaderBase() =
default;
114 bool use_internal_partition);
118 bool readGroupsFromFieldData(
IMesh*
mesh, vtkFieldData*,
int);
140 { VTK_FILE_EXT_VTU, &VtuMeshReaderBase::readMeshFromVtuFile },
149: m_trace_mng(trace_mng)
160bool VtuMeshReaderBase::
161readGroupsFromFieldData(
IMesh*
mesh, vtkFieldData* allFieldData,
int i)
163 const std::string group_name = allFieldData->GetArrayName(i);
164 vtkDataArray* iFieldData = allFieldData->GetArray(i);
168 if (iFieldData->GetDataType() != CURRENT_VTK_VERSION_LONG_TYPE)
171 auto nb_tuple = iFieldData->GetNumberOfTuples();
172 m_trace_mng->info() <<
"[readGroupsFromFieldData] iFieldData->GetNumberOfTuples=" << nb_tuple;
175 vtkLongArrayType* vtk_array = vtkLongArrayType::SafeDownCast(iFieldData);
179 IItemFamily* family = mesh->itemFamily(kind_type);
180 auto nb_item = nb_tuple - 1;
183 for (
Integer z = 0; z < nb_item; ++z)
184 unique_ids[z] = vtk_array->GetValue(z + 1);
188 family->itemsUniqueIdToLocalId(local_ids, unique_ids,
false);
193 for (
Integer index = 0; index < nb_item; ++index)
194 if (local_ids[index] != NULL_ITEM_LOCAL_ID)
195 ids.
add(local_ids[index]);
197 m_trace_mng->info() <<
"Create group family=" << family->name() <<
" name=" << group_name <<
" ids=" << ids.size();
198 family->createGroup(group_name, ids);
210 const String& dir_name,
bool use_internal_partition)
212 ARCANE_UNUSED(dir_name);
213 ARCANE_UNUSED(use_internal_partition);
214 bool itWasAnArcanProduction =
true;
215 IParallelMng* pm = mesh->parallelMng();
216 bool is_parallel = pm->isParallel();
217 Int32 sid = pm->commRank();
219 m_trace_mng->info() <<
"[readMeshFromVtuFile] Entering file_name=" << file_name;
220 vtkXMLUnstructuredGridReader* reader = vtkXMLUnstructuredGridReader::New();
221 std::string fname(file_name.toStdStringView());
222 reader->SetFileName(fname.c_str());
223 if (!reader->CanReadFile(fname.c_str())) {
224 m_trace_mng->info() <<
"Can not read file '" << file_name <<
"'";
227 reader->UpdateInformation();
229 vtkUnstructuredGrid* unstructuredGrid = reader->GetOutput();
232 auto nbOfCells = unstructuredGrid->GetNumberOfCells();
233 auto nbOfNodes = unstructuredGrid->GetNumberOfPoints();
238 m_trace_mng->info() <<
"[readMeshFromVtuFile] ## Now Fetching Nodes Unique IDs ##";
239 vtkPointData* allPointData = unstructuredGrid->GetPointData();
240 allPointData->Update();
241 vtkDataArray* dataNodeArray = allPointData->GetArray(
"NodesUniqueIDs");
242 vtkLongArrayType* nodesUidArray =
nullptr;
243 if (!dataNodeArray) {
244 m_trace_mng->info() <<
"[readMeshFromVtuFile] Could not be found, creating new one";
245 nodesUidArray = vtkLongArrayType::New();
246 for (
Integer uid = 0; uid < nbOfNodes; ++uid)
247 nodesUidArray->InsertNextValue(uid);
248 itWasAnArcanProduction =
false;
251 m_trace_mng->info() <<
"dataNodeArray->GetDataType()" << dataNodeArray->GetDataType() <<
" Long=" << CURRENT_VTK_VERSION_LONG_TYPE;
252 if (dataNodeArray->GetDataType() != CURRENT_VTK_VERSION_LONG_TYPE)
254 nodesUidArray = vtkLongArrayType::SafeDownCast(dataNodeArray);
256 m_trace_mng->info() <<
"[readMeshFromVtuFile] Fetched";
258 m_trace_mng->info() <<
"[readMeshFromVtuFile] ## Now Fetching Cells Unique IDs ##";
259 vtkCellData* allCellData = unstructuredGrid->GetCellData();
260 allCellData->Update();
261 vtkDataArray* dataCellArray = allCellData->GetArray(
"CellsUniqueIDs");
262 vtkLongArrayType* cellsUidArray =
nullptr;
263 if (!dataCellArray) {
264 m_trace_mng->info() <<
"[readMeshFromVtuFile] Could not be found, creating new one";
265 cellsUidArray = vtkLongArrayType::New();
266 for (
Integer uid = 0; uid < nbOfCells; ++uid)
267 cellsUidArray->InsertNextValue(uid);
268 itWasAnArcanProduction =
false;
271 if (dataCellArray->GetDataType() != CURRENT_VTK_VERSION_LONG_TYPE)
273 cellsUidArray = vtkLongArrayType::SafeDownCast(dataCellArray);
278 vtkDataArray* data_owner_array = allPointData->GetArray(
"NodesOwner");
279 vtkIntArray* nodes_owner_array =
nullptr;
280 if (data_owner_array) {
281 if (data_owner_array->GetDataType() == VTK_INT) {
282 nodes_owner_array = vtkIntArray::SafeDownCast(data_owner_array);
285 m_trace_mng->warning() <<
"NodesOwner array is present but has bad type (not Int32)";
288 m_trace_mng->info() <<
"[readMeshFromVtuFile] Fetched";
293 m_trace_mng->info() <<
"[readMeshFromVtuFile] nbOfCells=" << nbOfCells <<
", nbOfNodes=" << nbOfNodes;
295 HashTableMapT<Int64, Real3> nodes_coords(
Integer(nbOfNodes),
true);
296 vtkPoints* vtkAllPoints = unstructuredGrid->GetPoints();
298 mesh->toPrimaryMesh()->setDimension(vtkAllPoints->GetData()->GetNumberOfComponents());
299 for (vtkIdType i = 0; i < nbOfNodes; ++i) {
301 vtkAllPoints->GetPoint(i, xyz);
303 coords[i] = Real3(xyz[0], xyz[1], xyz[2]);
304 nodes_coords.nocheckAdd(nodesUidArray->GetValue(i), coords[i]);
308 HashTableMapT<Int64, Int32> nodes_owner_map(
Integer(nbOfNodes),
true);
309 if (nodes_owner_array) {
310 for (
Integer i = 0; i < nbOfNodes; ++i) {
311 nodes_owner_map.nocheckAdd(nodesUidArray->GetValue(i), nodes_owner_array->GetValue(i));
316 cells_filter.
resize(nbOfCells);
317 for (
Integer i = 0; i < nbOfCells; ++i)
321 vtkIdType mesh_nb_cell_node = 0;
322 for (
Integer j = 0, js = cells_filter.size(); j < js; ++j)
323 mesh_nb_cell_node += unstructuredGrid->GetCell(j)->GetNumberOfPoints();
324 m_trace_mng->info() <<
"Number of mesh_nb_cell_node = " << mesh_nb_cell_node;
334 for (
Integer i_cell = 0, s_cell = cells_filter.size(); i_cell < s_cell; ++i_cell) {
335 vtkIdType iVtkCell = i_cell;
336 vtkCell* cell = unstructuredGrid->GetCell(iVtkCell);
337 if (!cell->IsLinear())
338 throw NotImplementedException(A_FUNCINFO);
339 Integer arcItemType = IT_NullType;
340 switch (cell->GetCellType()) {
343 arcItemType = IT_Tetraedron4;
345 case (VTK_HEXAHEDRON):
346 arcItemType = IT_Hexaedron8;
349 arcItemType = IT_Pyramid5;
353 arcItemType = IT_Pentaedron6;
355 case (VTK_PENTAGONAL_PRISM):
356 arcItemType = IT_Heptaedron10;
358 case (VTK_HEXAGONAL_PRISM):
359 arcItemType = IT_Octaedron12;
363 if (mesh->dimension() == 2) {
364 arcItemType = IT_Quad4;
366 else if (mesh->dimension() == 3 && !mesh->meshKind().isMonoDimension()) {
367 arcItemType = IT_Cell3D_Quad4;
370 ARCANE_FATAL(
"VTK_QUAD is not supported in mono-dimension meshes of dimension {0}",
375 if (mesh->dimension() == 2) {
376 arcItemType = IT_Triangle3;
378 else if (mesh->dimension() == 3 && !mesh->meshKind().isMonoDimension()) {
379 arcItemType = IT_Cell3D_Triangle3;
382 ARCANE_FATAL(
"VTK_TRIANGLE is not supported in mono-dimension meshes of dimension {0}",
387 case (VTK_POLY_VERTEX):
388 switch (cell->GetNumberOfPoints()) {
390 arcItemType = IT_Tetraedron4;
393 arcItemType = IT_HemiHexa7;
396 arcItemType = IT_HemiHexa6;
399 arcItemType = IT_HemiHexa5;
402 throw NotImplementedException(A_FUNCINFO);
407 switch (cell->GetNumberOfPoints()) {
409 arcItemType = IT_Tetraedron4;
412 throw NotImplementedException(A_FUNCINFO);
416#ifndef NO_USER_WARNING
417#warning IT_Vertex returns 0 nodes in ItemTypeMng.cc@ 42:type->setInfos(this,IT_Vertex,0,0,0);
418#warning IT_Vertex vs IT_DualNode HACK
421 arcItemType = IT_DualNode;
424 switch (cell->GetNumberOfPoints()) {
426 if (mesh->dimension() == 2) {
427 arcItemType = IT_Triangle3;
429 else if (mesh->dimension() == 3 && !mesh->meshKind().isMonoDimension()) {
430 arcItemType = IT_Cell3D_Triangle3;
433 ARCANE_FATAL(
"VTK_TRIANGLE is not supported in mono-dimension meshes of dimension {0}",
438 if (mesh->dimension() == 2) {
439 arcItemType = IT_Quad4;
441 else if (mesh->dimension() == 3 && !mesh->meshKind().isMonoDimension()) {
442 arcItemType = IT_Cell3D_Quad4;
445 ARCANE_FATAL(
"VTK_QUAD is not supported in mono-dimension meshes of dimension {0}",
450 arcItemType = IT_Pentagon5;
453 arcItemType = IT_Hexagon6;
456 throw NotImplementedException(A_FUNCINFO);
461 case (VTK_POLY_LINE):
462 switch (cell->GetNumberOfPoints()) {
464 arcItemType = IT_CellLine2;
467 throw NotImplementedException(A_FUNCINFO);
471 case (VTK_EMPTY_CELL):
472 m_trace_mng->info() <<
"VTK_EMPTY_CELL n=" << cell->GetNumberOfPoints();
474 case (VTK_TRIANGLE_STRIP):
475 m_trace_mng->info() <<
"VTK_TRIANGLE_STRIP n=" << cell->GetNumberOfPoints();
478 m_trace_mng->info() <<
"VTK_VOXEL n=" << cell->GetNumberOfPoints();
479 switch (cell->GetNumberOfPoints()) {
481 arcItemType = IT_Triangle3;
484 arcItemType = IT_Quad4;
487 arcItemType = IT_Pentagon5;
490 arcItemType = IT_Hexagon6;
493 arcItemType = IT_HemiHexa7;
496 arcItemType = IT_Hexaedron8;
499 throw NotImplementedException(A_FUNCINFO);
504 throw NotImplementedException(A_FUNCINFO);
507 auto nNodes = cell->GetNumberOfPoints();
510 cells_infos[cells_infos_index++] = arcItemType;
513 Integer cell_indirect_id = cells_filter[i_cell];
514 cells_infos[cells_infos_index++] = cellsUidArray->GetValue(cell_indirect_id);
517 vtkIdList* nodeIds = cell->GetPointIds();
518 for (
int iNode = 0; iNode < nNodes; ++iNode) {
519 auto localId = nodeIds->GetId(iNode);
520 long uniqueUid = nodesUidArray->GetValue(localId);
522 cells_infos[cells_infos_index++] = uniqueUid;
529 m_trace_mng->info() <<
"[readMeshFromVtuFile] ## Mesh 3D ##";
530 PRIMARYMESH_CAST(mesh)->setDimension(3);
531 m_trace_mng->info() <<
"[readMeshFromVtuFile] ## Allocating " << cells_filter.size() <<
" cells ##";
532 PRIMARYMESH_CAST(mesh)->allocateCells(cells_filter.size(), cells_infos,
false);
536 m_trace_mng->info() <<
"[readMeshFromVtuFile] internalNodes.size()=" << internalNodes.size();
537 if (nodes_owner_array) {
538 m_trace_mng->info() <<
"Set nodes owners from vtu file";
539 for (
Integer i = 0, is = internalNodes.size(); i < is; ++i) {
540 ItemInternal* internal_node = internalNodes[i];
541 Int32 true_owner = nodes_owner_map[internal_node->uniqueId()];
542 internal_node->setOwner(true_owner, sid);
546 for (
Integer i = 0, is = internalNodes.size(); i < is; ++i)
547 internalNodes[i]->setOwner(sid, sid);
550 m_trace_mng->info() <<
"[readMeshFromVtuFile] internalCells.size()=" << internalCells.size();
551 for (
Integer i = 0, is = internalCells.size(); i < is; ++i)
552 internalCells[i]->setOwner(sid, sid);
557 m_trace_mng->info() <<
"[readMeshFromVtuFile] ## Ending with endAllocate ##";
558 PRIMARYMESH_CAST(mesh)->endAllocate();
561 mesh->utilities()->changeOwnersFromCells();
564 m_trace_mng->info() <<
"\n[readMeshFromVtuFile] ## Now dealing with ghost's layer ##";
565 m_trace_mng->info() <<
"[readMeshFromVtuFile] mesh.nbNode=" << mesh->nbNode() <<
" mesh.nbCell=" << mesh->nbCell();
570 m_trace_mng->info() <<
"[readMeshFromVtuFile] ## Now Fetching Other Fields ##";
571 vtkFieldData* allFieldData = unstructuredGrid->GetFieldData();
572 int nbrOfArrays = allFieldData->GetNumberOfArrays();
573 m_trace_mng->info() <<
"[readMeshFromVtuFile] nbrOfArrays = " << nbrOfArrays;
574 for (
int i = 0; i < nbrOfArrays; ++i) {
575 if (itWasAnArcanProduction ==
false)
577 m_trace_mng->info() <<
"[readMeshFromVtuFile] Focussing on \"" << allFieldData->GetArrayName(i) <<
"\" (i="
579 if (readGroupsFromFieldData(mesh, allFieldData, i) !=
true)
586 m_trace_mng->info() <<
"[readMeshFromVtuFile] ## Now insert coords ##";
590 nodes_coord_var[inode] = nodes_coords[inode->uniqueId()];
596 mesh->synchronizeGroupsAndVariables();
602 m_trace_mng->info() <<
"[readMeshFromVtuFile] RTOk";
606 nodesUidArray->Delete();
608 cellsUidArray->Delete();
635 bool use_internal_partition)
override;
651 for (vtkFileExtIdx = 0; vtkFileExtReader[vtkFileExtIdx].ext !=
nullptr; ++vtkFileExtIdx) {
653 if (str == vtkFileExtReader[vtkFileExtIdx].ext) {
672 bool use_internal_partition)
674 ARCANE_UNUSED(mesh_node);
675 info() <<
"[readMeshFromFile] Forwarding to vtkFileExtReader[" << vtkFileExtIdx <<
"].reader";
678 return vtu_mesh_reader.readMeshFromVtuFile(
mesh, file_name, dir_name, use_internal_partition);
696class VtuCaseMeshReader
709 , m_read_info(read_info)
717 ARCANE_UNUSED(build_info);
723 String fname = m_read_info.fileName();
724 m_trace_mng->info() <<
"Vtu Reader (ICaseMeshReader) file_name=" << fname;
725 bool ret = vtu_service.readMeshFromVtuFile(pm, fname, m_read_info.directoryName(),
726 m_read_info.isParallelRead());
749 if (read_info.format() ==
"vtu")
#define ARCANE_FATAL(...)
Macro throwing a FatalErrorException.
#define ARCANE_REGISTER_SUB_DOMAIN_FACTORY(aclass, ainterface, aname)
Registers a factory service for the class aclass.
#define ARCANE_SERVICE_INTERFACE(ainterface)
Macro to declare an interface when registering a service.
AbstractService(const ServiceBuildInfo &)
Constructor from a ServiceBuildInfo.
void resize(Int64 s)
Changes the number of elements in the array to s.
void add(ConstReferenceType val)
Adds element val to the end of the array.
Necessary information for reading a mesh file.
Interface for the mesh reading service from the dataset.
virtual ITraceMng * traceMng()=0
Associated message manager.
Interface of a mesh creation/reading service.
Interface of the service managing the reading of a mesh.
eReturnType
Types of return codes for a read or write operation.
@ RTError
Error during the operation.
@ RTOk
Operation successfully performed.
Interface of the subdomain manager.
Parameters necessary for building a mesh.
Reference to an instance.
ISubDomain * subDomain() const
Access to the associated ISubDomain.
Structure containing the information to create a service.
Service creation properties.
Unicode character string.
TraceMessage info() const
Flow for an information message.
ITraceMng * traceMng() const
Trace manager.
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.
Reader for mesh files in VTK format.
bool allowExtension(const String &str) override
Checks if the service supports files with the extension str.
eReturnType readMeshFromFile(IPrimaryMesh *mesh, const XmlNode &mesh_node, const String &file_name, const String &dir_name, bool use_internal_partition) override
Reads a mesh from a file.
#define ARCANE_REGISTER_SERVICE(aclass, a_service_property,...)
Macro for registering a service.
MeshVariableScalarRefT< Node, Real3 > VariableNodeReal3
Coordinate type quantity at node.
-- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature --
UniqueArray< Int64 > Int64UniqueArray
Dynamic 1D array of 64-bit integers.
Int32 Integer
Type representing an integer.
UniqueArray< Real3 > Real3UniqueArray
Dynamic 1D array of rank 3 vectors.
ConstArrayView< ItemInternal * > ItemInternalList
Type of the internal list of entities.
@ ST_SubDomain
The service is used at the subdomain level.
UniqueArray< Int32 > Int32UniqueArray
Dynamic 1D array of 32-bit integers.
eItemKind
Mesh entity type.
@ IK_Node
Node mesh entity.
@ IK_Cell
Cell mesh entity.
auto makeRef(InstanceType *t) -> Ref< InstanceType >
Creates a reference on a pointer.
UniqueArray< Integer > IntegerUniqueArray
Dynamic 1D array of integers.
std::int32_t Int32
Signed integer type of 32 bits.