Arcane  4.1.12.0
Developer documentation
Loading...
Searching...
No Matches
VtuMeshReader.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2026 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* VtuMeshReader.cc (C) 2000-2025 */
9/* */
10/* Reading/Writing a file in VtuMeshReader format. */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
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"
27
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"
53
54#include <vtkXMLUnstructuredGridReader.h>
55#include <vtkUnstructuredGrid.h>
56#include <vtkCell.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>
63
64/**********************************************************************
65 * Convention data types/file type and particular file extension parity *
66 **********************************************************************/
67#define VTK_FILE_EXT_VTI "vti" // Serial structured vtkImageData
68#define VTK_FILE_EXT_VTP "vtp" // Serial UNstructured vtkPolyData
69#define VTK_FILE_EXT_VTR "vtr" // Serial structured vtkRectilinearGrid
70#define VTK_FILE_EXT_VTS "vts" // Serial structured vtkStructuredGrid
71#define VTK_FILE_EXT_VTU "vtu" // Serial UNstructured vtkUnstructuredGrid
72#define VTK_FILE_EXT_PVTI "pvti" // Parallel structured vtkImageData
73#define VTK_FILE_EXT_PVTP "pvtp" // Parallel UNstructured vtkPolyData
74#define VTK_FILE_EXT_PVTR "pvtr" // Parallel structured vtkRectilinearGrid
75#define VTK_FILE_EXT_PVTS "pvts" // Parallel structured vtkStructuredGrid
76#define VTK_FILE_EXT_PVTU "pvtu" // Parallel UNstructured vtkUnstructuredGrid
77
78#if VTK_MAJOR_VERSION >= 9
79#define CURRENT_VTK_VERSION_LONG_TYPE VTK_LONG_LONG
80#include <vtkLongLongArray.h>
81using vtkLongArrayType = vtkLongLongArray;
82#else
83#define CURRENT_VTK_VERSION_LONG_TYPE VTK_LONG
84using vtkLongArrayType = vtkLongArray;
85#endif
86
87/*---------------------------------------------------------------------------*/
88/*---------------------------------------------------------------------------*/
89
90namespace Arcane
91{
92
93/*---------------------------------------------------------------------------*/
94/*---------------------------------------------------------------------------*/
95
99class VtuMeshReaderBase
100{
101 public:
102
103 explicit VtuMeshReaderBase(ITraceMng* trace_mng);
104 virtual ~VtuMeshReaderBase() = default;
105
106 public:
107
108 virtual void build()
109 {
110 }
111
112 IMeshReader::eReturnType readMeshFromVtuFile(IMesh* mesh,
113 const String& file_name, const String& dir_name,
114 bool use_internal_partition);
115
116 // ISubDomain* subDomain() { return m_trace_mng; }
117
118 bool readGroupsFromFieldData(IMesh* mesh, vtkFieldData*, int);
119
120 private:
121
122 int vtkFileExtIdx;
123 ITraceMng* m_trace_mng;
124 // IMeshReader::eReturnType writeMeshToLemFile(IMesh* mesh, const String& filename,const String& dir_name);
125};
126
127/*---------------------------------------------------------------------------*/
128/*---------------------------------------------------------------------------*/
129
130/****************************************************************************
131 * Assignation des data types/file-type and particular file extension readers *
132\****************************************************************************/
133typedef struct
134{
135 char* ext;
136 IMeshReader::eReturnType (VtuMeshReaderBase::*reader)(IMesh*, const String&, const String&, bool);
138
139vtkExtReader vtkFileExtReader[] = {
140 { VTK_FILE_EXT_VTU, &VtuMeshReaderBase::readMeshFromVtuFile },
141 { nullptr, nullptr }
142};
143
144/*---------------------------------------------------------------------------*/
145/*---------------------------------------------------------------------------*/
146
147VtuMeshReaderBase::
148VtuMeshReaderBase(ITraceMng* trace_mng)
149: m_trace_mng(trace_mng)
150{
151 vtkFileExtIdx = 0;
152}
153
154/*---------------------------------------------------------------------------*/
155/*---------------------------------------------------------------------------*/
156
157/*****************************************************************
158 * readGroupsFromFieldData
159 *****************************************************************/
160bool VtuMeshReaderBase::
161readGroupsFromFieldData(IMesh* mesh, vtkFieldData* allFieldData, int i)
162{
163 const std::string group_name = allFieldData->GetArrayName(i);
164 vtkDataArray* iFieldData = allFieldData->GetArray(i);
165 if (!iFieldData)
166 return false;
167
168 if (iFieldData->GetDataType() != CURRENT_VTK_VERSION_LONG_TYPE)
169 return false;
170
171 auto nb_tuple = iFieldData->GetNumberOfTuples();
172 m_trace_mng->info() << "[readGroupsFromFieldData] iFieldData->GetNumberOfTuples=" << nb_tuple;
173 if (nb_tuple == 0)
174 return false;
175 vtkLongArrayType* vtk_array = vtkLongArrayType::SafeDownCast(iFieldData);
176
177 // The first element of the array contains its type
178 eItemKind kind_type = (eItemKind)(vtk_array->GetValue(0));
179 IItemFamily* family = mesh->itemFamily(kind_type);
180 auto nb_item = nb_tuple - 1;
181 Int64UniqueArray unique_ids(nb_item);
182 // The following elements contain the uniqueId() of the group entities.
183 for (Integer z = 0; z < nb_item; ++z)
184 unique_ids[z] = vtk_array->GetValue(z + 1);
185
186 // Retrieves the corresponding localId().
187 Int32UniqueArray local_ids(unique_ids.size());
188 family->itemsUniqueIdToLocalId(local_ids, unique_ids, false);
189
190 // Not all entities are necessarily in the current mesh and
191 // they must therefore be filtered.
193 for (Integer index = 0; index < nb_item; ++index)
194 if (local_ids[index] != NULL_ITEM_LOCAL_ID)
195 ids.add(local_ids[index]);
196
197 m_trace_mng->info() << "Create group family=" << family->name() << " name=" << group_name << " ids=" << ids.size();
198 family->createGroup(group_name, ids);
199 return true;
200}
201
202/*---------------------------------------------------------------------------*/
203/*---------------------------------------------------------------------------*/
204
205/*****************************************************************
206 * readMeshFromVtuFile
207 *****************************************************************/
208IMeshReader::eReturnType VtuMeshReaderBase::
209readMeshFromVtuFile(IMesh* mesh, const String& file_name,
210 const String& dir_name, bool use_internal_partition)
211{
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();
218
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 << "'";
226 }
227 reader->UpdateInformation();
228 reader->Update(); // Force reading
229 vtkUnstructuredGrid* unstructuredGrid = reader->GetOutput();
230 // With VTK 7, no longer need the Update.
231 //unstructuredGrid->Update();// The actual file reading only takes place after calling Update().
232 auto nbOfCells = unstructuredGrid->GetNumberOfCells();
233 auto nbOfNodes = unstructuredGrid->GetNumberOfPoints();
234
235 /*******************************
236 *Fetching Nodes UID & Cells UID *
237 *******************************/
238 m_trace_mng->info() << "[readMeshFromVtuFile] ## Now Fetching Nodes Unique IDs ##";
239 vtkPointData* allPointData = unstructuredGrid->GetPointData(); // Data associated to Points
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;
249 }
250 else {
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);
255 }
256 m_trace_mng->info() << "[readMeshFromVtuFile] Fetched";
257
258 m_trace_mng->info() << "[readMeshFromVtuFile] ## Now Fetching Cells Unique IDs ##";
259 vtkCellData* allCellData = unstructuredGrid->GetCellData(); // Data associated to Points
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;
269 }
270 else {
271 if (dataCellArray->GetDataType() != CURRENT_VTK_VERSION_LONG_TYPE)
273 cellsUidArray = vtkLongArrayType::SafeDownCast(dataCellArray);
274 }
275
276 // Table containing the owner numbers of the nodes.
277 // This table is optional and is used by the partitioner
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);
283 }
284 else {
285 m_trace_mng->warning() << "NodesOwner array is present but has bad type (not Int32)";
286 }
287 }
288 m_trace_mng->info() << "[readMeshFromVtuFile] Fetched";
289
290 /************************
291 * Fetch own nodes coords *
292 ************************/
293 m_trace_mng->info() << "[readMeshFromVtuFile] nbOfCells=" << nbOfCells << ", nbOfNodes=" << nbOfNodes;
294 Real3UniqueArray coords(nbOfNodes);
295 HashTableMapT<Int64, Real3> nodes_coords(Integer(nbOfNodes), true);
296 vtkPoints* vtkAllPoints = unstructuredGrid->GetPoints();
297 // Get mesh dimension
298 mesh->toPrimaryMesh()->setDimension(vtkAllPoints->GetData()->GetNumberOfComponents());
299 for (vtkIdType i = 0; i < nbOfNodes; ++i) {
300 double xyz[3];
301 vtkAllPoints->GetPoint(i, xyz);
302 //m_trace_mng->info() << "x=" << xyz[0] << " y=" << xyz[1]<< " z=" << xyz[2] << ", nodUid=" << nodesUidArray->GetValue(i);
303 coords[i] = Real3(xyz[0], xyz[1], xyz[2]);
304 nodes_coords.nocheckAdd(nodesUidArray->GetValue(i), coords[i]);
305 }
306
307 // Create hash table for nodes owner.
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));
312 }
313 }
314
315 IntegerUniqueArray cells_filter;
316 cells_filter.resize(nbOfCells);
317 for (Integer i = 0; i < nbOfCells; ++i)
318 cells_filter[i] = i;
319
320 // Calculate the number of cells/nodes
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;
325
326 // Table containing the info for the cells (see IMesh::allocateMesh())
327 Int64UniqueArray cells_infos(mesh_nb_cell_node + cells_filter.size() * 2);
328 Integer cells_infos_index = 0;
329
330 /******************
331 * Now insert CELLS *
332 ******************/
333 // For all cells that are discovered
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()) {
341 // Linear cells
342 case (VTK_TETRA):
343 arcItemType = IT_Tetraedron4;
344 break;
345 case (VTK_HEXAHEDRON):
346 arcItemType = IT_Hexaedron8;
347 break;
348 case (VTK_PYRAMID):
349 arcItemType = IT_Pyramid5;
350 break;
351 /* Others not yet implemented */
352 case (VTK_WEDGE):
353 arcItemType = IT_Pentaedron6;
354 break;
355 case (VTK_PENTAGONAL_PRISM):
356 arcItemType = IT_Heptaedron10;
357 break;
358 case (VTK_HEXAGONAL_PRISM):
359 arcItemType = IT_Octaedron12;
360 break;
361 /* 2D */
362 case (VTK_QUAD):
363 if (mesh->dimension() == 2) {
364 arcItemType = IT_Quad4;
365 }
366 else if (mesh->dimension() == 3 && !mesh->meshKind().isMonoDimension()) {
367 arcItemType = IT_Cell3D_Quad4;
368 }
369 else {
370 ARCANE_FATAL("VTK_QUAD is not supported in mono-dimension meshes of dimension {0}",
371 mesh->dimension());
372 }
373 break;
374 case (VTK_TRIANGLE):
375 if (mesh->dimension() == 2) {
376 arcItemType = IT_Triangle3;
377 }
378 else if (mesh->dimension() == 3 && !mesh->meshKind().isMonoDimension()) {
379 arcItemType = IT_Cell3D_Triangle3;
380 }
381 else {
382 ARCANE_FATAL("VTK_TRIANGLE is not supported in mono-dimension meshes of dimension {0}",
383 mesh->dimension());
384 }
385 break;
386 // Case for parsing poly vertex
387 case (VTK_POLY_VERTEX):
388 switch (cell->GetNumberOfPoints()) {
389 case (4):
390 arcItemType = IT_Tetraedron4;
391 break;
392 case (7):
393 arcItemType = IT_HemiHexa7;
394 break;
395 case (6):
396 arcItemType = IT_HemiHexa6;
397 break;
398 case (5):
399 arcItemType = IT_HemiHexa5;
400 break;
401 default:
402 throw NotImplementedException(A_FUNCINFO);
403 ;
404 }
405 break;
406 case (VTK_PIXEL):
407 switch (cell->GetNumberOfPoints()) {
408 case (4):
409 arcItemType = IT_Tetraedron4;
410 break;
411 default:
412 throw NotImplementedException(A_FUNCINFO);
413 ;
414 }
415 break;
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
419#endif /* NO_USER_WARNING */
420 case (VTK_VERTEX):
421 arcItemType = IT_DualNode;
422 break;
423 case (VTK_POLYGON):
424 switch (cell->GetNumberOfPoints()) {
425 case (3):
426 if (mesh->dimension() == 2) {
427 arcItemType = IT_Triangle3;
428 }
429 else if (mesh->dimension() == 3 && !mesh->meshKind().isMonoDimension()) {
430 arcItemType = IT_Cell3D_Triangle3;
431 }
432 else {
433 ARCANE_FATAL("VTK_TRIANGLE is not supported in mono-dimension meshes of dimension {0}",
434 mesh->dimension());
435 }
436 break;
437 case (4):
438 if (mesh->dimension() == 2) {
439 arcItemType = IT_Quad4;
440 }
441 else if (mesh->dimension() == 3 && !mesh->meshKind().isMonoDimension()) {
442 arcItemType = IT_Cell3D_Quad4;
443 }
444 else {
445 ARCANE_FATAL("VTK_QUAD is not supported in mono-dimension meshes of dimension {0}",
446 mesh->dimension());
447 }
448 break;
449 case (5):
450 arcItemType = IT_Pentagon5;
451 break;
452 case (6):
453 arcItemType = IT_Hexagon6;
454 break;
455 default:
456 throw NotImplementedException(A_FUNCINFO);
457 ;
458 }
459 break;
460 case (VTK_LINE):
461 case (VTK_POLY_LINE):
462 switch (cell->GetNumberOfPoints()) {
463 case (2):
464 arcItemType = IT_CellLine2;
465 break;
466 default:
467 throw NotImplementedException(A_FUNCINFO);
468 ;
469 }
470 break;
471 case (VTK_EMPTY_CELL):
472 m_trace_mng->info() << "VTK_EMPTY_CELL n=" << cell->GetNumberOfPoints();
473 break;
474 case (VTK_TRIANGLE_STRIP):
475 m_trace_mng->info() << "VTK_TRIANGLE_STRIP n=" << cell->GetNumberOfPoints();
476 break;
477 case (VTK_VOXEL):
478 m_trace_mng->info() << "VTK_VOXEL n=" << cell->GetNumberOfPoints();
479 switch (cell->GetNumberOfPoints()) {
480 case (3):
481 arcItemType = IT_Triangle3;
482 break;
483 case (4):
484 arcItemType = IT_Quad4;
485 break;
486 case (5):
487 arcItemType = IT_Pentagon5;
488 break;
489 case (6):
490 arcItemType = IT_Hexagon6;
491 break;
492 case (7):
493 arcItemType = IT_HemiHexa7;
494 break;
495 case (8):
496 arcItemType = IT_Hexaedron8;
497 break;
498 default:
499 throw NotImplementedException(A_FUNCINFO);
500 ;
501 }
502 break;
503 default:
504 throw NotImplementedException(A_FUNCINFO);
505 ;
506 }
507 auto nNodes = cell->GetNumberOfPoints(); // Return the number of points in the cell
508
509 // First is cell's TYPE Stores the cell type
510 cells_infos[cells_infos_index++] = arcItemType;
511
512 // Then comes its UniqueID Stores the cell unique ID
513 Integer cell_indirect_id = cells_filter[i_cell];
514 cells_infos[cells_infos_index++] = cellsUidArray->GetValue(cell_indirect_id);
515
516 // And finally the Nodes' unique IDs
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);
521 //info() << "working on localId=" << localId << ", uniqueUid=" << uniqueUid;
522 cells_infos[cells_infos_index++] = uniqueUid;
523 }
524 }
525
526 /********************************
527 * Setting Dimension & Allocating *
528 ********************************/
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);
533
534 // Positions the node owners from the node groups
535 ItemInternalList internalNodes(mesh->itemsInternal(IK_Node));
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);
543 }
544 }
545 else
546 for (Integer i = 0, is = internalNodes.size(); i < is; ++i)
547 internalNodes[i]->setOwner(sid, sid);
548
549 ItemInternalList internalCells(mesh->itemsInternal(IK_Cell));
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);
553
554 /********************************************
555 * Now finishing & preparing for ghost layout *
556 ********************************************/
557 m_trace_mng->info() << "[readMeshFromVtuFile] ## Ending with endAllocate ##";
558 PRIMARYMESH_CAST(mesh)->endAllocate();
559 if (is_parallel) {
560 // mesh->setOwnersFromCells();
561 mesh->utilities()->changeOwnersFromCells();
562 }
563
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();
566
567 /***********************
568 * Fetching Other Groups *
569 ***********************/
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)
576 continue;
577 m_trace_mng->info() << "[readMeshFromVtuFile] Focussing on \"" << allFieldData->GetArrayName(i) << "\" (i="
578 << i << ")";
579 if (readGroupsFromFieldData(mesh, allFieldData, i) != true)
581 }
582
583 /*******************
584 * Now insert coords *
585 *******************/
586 m_trace_mng->info() << "[readMeshFromVtuFile] ## Now insert coords ##";
587 // Fills the variable containing the node coordinates
588 VariableNodeReal3& nodes_coord_var(PRIMARYMESH_CAST(mesh)->nodesCoordinates());
589 ENUMERATE_NODE (inode, mesh->ownNodes()) {
590 nodes_coord_var[inode] = nodes_coords[inode->uniqueId()];
591 }
592
593 /****************************************
594 * Synchronizing groups/variables & nodes *
595 ****************************************/
596 mesh->synchronizeGroupsAndVariables();
597
598 /**************
599 * Finishing up *
600 **************/
601 reader->Delete();
602 m_trace_mng->info() << "[readMeshFromVtuFile] RTOk";
603
604 // TODO: look into automatic destruction
605 if (!dataNodeArray)
606 nodesUidArray->Delete();
607 if (!dataCellArray)
608 cellsUidArray->Delete();
609
610 return IMeshReader::RTOk;
611}
612
613/*---------------------------------------------------------------------------*/
614/*---------------------------------------------------------------------------*/
615
616class VtuMeshReader
617: public AbstractService
618, public IMeshReader
619{
620 public:
621
622 explicit VtuMeshReader(const ServiceBuildInfo& sbi)
623 : AbstractService(sbi)
624 , vtkFileExtIdx(-1)
625 , m_sub_domain(sbi.subDomain())
626 {
627 }
628
629 public:
630
631 bool allowExtension(const String& str) override;
632
634 const String& file_name, const String& dir_name,
635 bool use_internal_partition) override;
636
637 private:
638
639 int vtkFileExtIdx;
640 ISubDomain* m_sub_domain;
641 // eReturnType writeMeshToLemFile(IMesh* mesh, const String& filename,const String& dir_name);
642};
643
644/*****************************************************************
645 * allowExtension
646 *****************************************************************/
648allowExtension(const String& str)
649{
650 //m_trace_mng->info() << "[allowExtension] Checking for file extension...";
651 for (vtkFileExtIdx = 0; vtkFileExtReader[vtkFileExtIdx].ext != nullptr; ++vtkFileExtIdx) {
652 //m_trace_mng->info() << "Testing for '" << vtkFileExtReader[vtkFileExtIdx].ext << "'...";
653 if (str == vtkFileExtReader[vtkFileExtIdx].ext) {
654 return true;
655 }
656 }
657 //m_trace_mng->info()<<"Miss for all!";
658 // Sets our index in place for further service management
659 vtkFileExtIdx = 0;
660 return false;
661}
662
663/*---------------------------------------------------------------------------*/
664/*---------------------------------------------------------------------------*/
665
666/*****************************************************************
667 * readMeshFromFile switch
668 *****************************************************************/
670readMeshFromFile(IPrimaryMesh* mesh, const XmlNode& mesh_node,
671 const String& file_name, const String& dir_name,
672 bool use_internal_partition)
673{
674 ARCANE_UNUSED(mesh_node);
675 info() << "[readMeshFromFile] Forwarding to vtkFileExtReader[" << vtkFileExtIdx << "].reader";
676 VtuMeshReaderBase vtu_mesh_reader{ traceMng() };
677 // return (vtu_mesh_reader.*vtkFileExtReader[vtkFileExtIdx].reader)(mesh, file_name, dir_name, use_internal_partition);
678 return vtu_mesh_reader.readMeshFromVtuFile(mesh, file_name, dir_name, use_internal_partition);
679}
680
681/*---------------------------------------------------------------------------*/
682/*---------------------------------------------------------------------------*/
683
685
687 ServiceProperty("VtuNewMeshReader", ST_SubDomain),
689
690/*---------------------------------------------------------------------------*/
691/*---------------------------------------------------------------------------*/
692
693/*---------------------------------------------------------------------------*/
694/*---------------------------------------------------------------------------*/
695
696class VtuCaseMeshReader
697: public AbstractService
698, public ICaseMeshReader
699{
700 public:
701
702 class Builder
703 : public IMeshBuilder
704 {
705 public:
706
707 explicit Builder(ITraceMng* tm, const CaseMeshReaderReadInfo& read_info)
708 : m_trace_mng(tm)
709 , m_read_info(read_info)
710 {
711 }
712
713 public:
714
715 void fillMeshBuildInfo(MeshBuildInfo& build_info) override
716 {
717 ARCANE_UNUSED(build_info);
718 }
719
721 {
722 VtuMeshReaderBase vtu_service(pm->traceMng());
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());
727 if (ret)
728 ARCANE_FATAL("Can not read VTK File");
729 }
730
731 private:
732
733 ITraceMng* m_trace_mng;
734 CaseMeshReaderReadInfo m_read_info;
735 };
736
737 public:
738
739 explicit VtuCaseMeshReader(const ServiceBuildInfo& sbi)
740 : AbstractService(sbi)
741 {
742 }
743
744 public:
745
747 {
748 IMeshBuilder* builder = nullptr;
749 if (read_info.format() == "vtu")
750 builder = new Builder(traceMng(), read_info);
751 return makeRef(builder);
752 }
753};
754
755/*---------------------------------------------------------------------------*/
756/*---------------------------------------------------------------------------*/
757
758ARCANE_REGISTER_SERVICE(VtuCaseMeshReader,
759 ServiceProperty("VtuCaseMeshReader", ST_SubDomain),
760 ARCANE_SERVICE_INTERFACE(ICaseMeshReader));
761
762/*---------------------------------------------------------------------------*/
763/*---------------------------------------------------------------------------*/
764
765} // End namespace Arcane
766
767/*---------------------------------------------------------------------------*/
768/*---------------------------------------------------------------------------*/
#define ARCANE_FATAL(...)
Macro throwing a FatalErrorException.
Types and macros for iterating over mesh entities.
#define ENUMERATE_NODE(name, group)
Generic enumerator for a node group.
#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.
Base class of 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.
Definition IMeshReader.h:33
eReturnType
Types of return codes for a read or write operation.
Definition IMeshReader.h:38
@ RTError
Error during the operation.
Definition IMeshReader.h:40
@ RTOk
Operation successfully performed.
Definition IMeshReader.h:39
Interface of the subdomain manager.
Definition ISubDomain.h:75
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.
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.
Node of a DOM tree.
Definition XmlNode.h:51
#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.
Definition UtilsTypes.h:339
Int32 Integer
Type representing an integer.
UniqueArray< Real3 > Real3UniqueArray
Dynamic 1D array of rank 3 vectors.
Definition UtilsTypes.h:363
ConstArrayView< ItemInternal * > ItemInternalList
Type of the internal list of entities.
Definition ItemTypes.h:466
@ ST_SubDomain
The service is used at the subdomain level.
UniqueArray< Int32 > Int32UniqueArray
Dynamic 1D array of 32-bit integers.
Definition UtilsTypes.h:341
eItemKind
Mesh entity type.
@ IK_Node
Node mesh entity.
@ IK_Cell
Cell mesh entity.
auto makeRef(InstanceType *t) -> Ref< InstanceType >
Creates a reference on a pointer.
UniqueArray< Integer > IntegerUniqueArray
Dynamic 1D array of integers.
Definition UtilsTypes.h:347
std::int32_t Int32
Signed integer type of 32 bits.