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/AbstractService.h"
29#include "arcane/FactoryService.h"
30#include "arcane/IMainFactory.h"
31#include "arcane/IMeshReader.h"
32#include "arcane/ISubDomain.h"
33#include "arcane/IMesh.h"
34#include "arcane/IMeshSubMeshTransition.h"
35#include "arcane/IItemFamily.h"
36#include "arcane/Item.h"
37#include "arcane/ItemEnumerator.h"
38#include "arcane/VariableTypes.h"
39#include "arcane/IVariableAccessor.h"
40#include "arcane/IParallelMng.h"
41#include "arcane/IIOMng.h"
42#include "arcane/IXmlDocumentHolder.h"
43#include "arcane/XmlNodeList.h"
44#include "arcane/XmlNode.h"
45#include "arcane/IMeshUtilities.h"
46#include "arcane/IMeshWriter.h"
47#include "arcane/BasicService.h"
48#include "arcane/ItemPrinter.h"
49#include "arcane/ServiceBuilder.h"
51#include <vtkXMLUnstructuredGridReader.h>
52#include <vtkUnstructuredGrid.h>
54#include <vtkDataObject.h>
55#include <vtkDataSet.h>
56#include <vtkCellData.h>
57#include <vtkPointData.h>
58#include <vtkLongArray.h>
59#include <vtkIntArray.h>
65#define VTK_FILE_EXT_VTI "vti"
66#define VTK_FILE_EXT_VTP "vtp"
67#define VTK_FILE_EXT_VTR "vtr"
68#define VTK_FILE_EXT_VTS "vts"
69#define VTK_FILE_EXT_VTU "vtu"
70#define VTK_FILE_EXT_PVTI "pvti"
71#define VTK_FILE_EXT_PVTP "pvtp"
72#define VTK_FILE_EXT_PVTR "pvtr"
73#define VTK_FILE_EXT_PVTS "pvts"
74#define VTK_FILE_EXT_PVTU "pvtu"
102 bool allowExtension(
const String& str);
109 ISubDomain* subDomain() {
return m_sub_domain; }
111 bool readGroupsFromFieldData(IMesh *mesh, vtkFieldData*,
int);
116 ISubDomain* m_sub_domain;
117 eReturnType writeMeshToLemFile(IMesh* mesh,
const String& filename,
const String& dir_name);
137 {VTK_FILE_EXT_VTU, &VtuMeshReader::readMeshFromVtuFile},
146VtuMeshReader(
const ServiceBuildInfo& sbi)
147: AbstractService(sbi)
148, m_sub_domain(sbi.subDomain())
163 for(vtkFileExtIdx=0;vtkFileExtReader[vtkFileExtIdx].ext!=
NULL;++vtkFileExtIdx){
165 if (str == vtkFileExtReader[vtkFileExtIdx].ext){
186 info() <<
"[readMeshFromFile] Forwarding to vtkFileExtReader[" << vtkFileExtIdx <<
"].reader";
204 pfatal() <<
"Mesh writer service selected not available";
218readGroupsFromFieldData(IMesh* mesh,vtkFieldData* allFieldData,
int i)
220 const std::string group_name = allFieldData->GetArrayName(i);
221 vtkDataArray* iFieldData = allFieldData->GetArray(i);
225 if (iFieldData->GetDataType() != VTK_LONG)
228 Integer nb_tuple = iFieldData->GetNumberOfTuples();
229 info() <<
"[readGroupsFromFieldData] iFieldData->GetNumberOfTuples="<< nb_tuple;
232 vtkLongArray* vtk_array = vtkLongArray::SafeDownCast(iFieldData);
236 IItemFamily* family = mesh->itemFamily(kind_type);
237 Integer nb_item = nb_tuple - 1;
240 for( Integer z=0; z<nb_item; ++z )
241 unique_ids[z] = vtk_array->GetValue(z+1);
245 family->itemsUniqueIdToLocalId(local_ids,unique_ids,
false);
250 for( Integer i=0; i<nb_item; ++i )
251 if (local_ids[i]!=NULL_ITEM_LOCAL_ID)
252 ids.add(local_ids[i]);
254 info() <<
"Create group family=" << family->name() <<
" name=" << group_name <<
" ids=" << ids.size();
255 family->createGroup(group_name,ids);
266readMeshFromVtuFile(IMesh* mesh,
const XmlNode& mesh_node,
267 const String& filename,
const String& dir_name,
268 bool use_internal_partition)
270 bool itWasAnArcanProduction=
true;
271 IParallelMng* pm = mesh->parallelMng();
275 info() <<
"[readMeshFromVtuFile] Entering";
276 vtkXMLUnstructuredGridReader *reader = vtkXMLUnstructuredGridReader::New();
277 std::string fname = filename.localstr();
278 reader->SetFileName(fname.c_str());
279 if (!reader->CanReadFile(fname.c_str()))
281 reader->UpdateInformation();
283 vtkUnstructuredGrid *unstructuredGrid = reader->GetOutput();
286 int nbOfCells = unstructuredGrid->GetNumberOfCells();
287 int nbOfNodes = unstructuredGrid->GetNumberOfPoints();
293 info() <<
"[readMeshFromVtuFile] ## Now Fetching Nodes Unique IDs ##";
294 vtkPointData* allPointData=unstructuredGrid->GetPointData();
295 allPointData->Update();
296 vtkDataArray* dataNodeArray = allPointData->GetArray(
"NodesUniqueIDs");
297 vtkLongArray *nodesUidArray =
nullptr;
299 info() <<
"[readMeshFromVtuFile] Could not be found, creating new one";
300 nodesUidArray = vtkLongArray::New();
301 for(Integer uid=0; uid<nbOfNodes; ++uid)
302 nodesUidArray->InsertNextValue(uid);
303 itWasAnArcanProduction=
false;
306 if (dataNodeArray->GetDataType() != VTK_LONG)
308 nodesUidArray = vtkLongArray::SafeDownCast(dataNodeArray);
310 info() <<
"[readMeshFromVtuFile] Fetched";
312 info() <<
"[readMeshFromVtuFile] ## Now Fetching Cells Unique IDs ##";
313 vtkCellData* allCellData=unstructuredGrid->GetCellData();
314 allCellData->Update();
315 vtkDataArray* dataCellArray = allCellData->GetArray(
"CellsUniqueIDs");
316 vtkLongArray *cellsUidArray =
nullptr;
318 info() <<
"[readMeshFromVtuFile] Could not be found, creating new one";
319 cellsUidArray = vtkLongArray::New();
320 for(Integer uid=0; uid<nbOfCells; ++uid)
321 cellsUidArray->InsertNextValue(uid);
322 itWasAnArcanProduction=
false;
325 if (dataCellArray->GetDataType() != VTK_LONG)
327 cellsUidArray = vtkLongArray::SafeDownCast(dataCellArray);
332 vtkDataArray* data_owner_array = allPointData->GetArray(
"NodesOwner");
333 vtkIntArray* nodes_owner_array = 0;
334 if (data_owner_array){
335 if (data_owner_array->GetDataType() == VTK_INT){
336 nodes_owner_array = vtkIntArray::SafeDownCast(data_owner_array);
339 warning() <<
"NodesOwner array is present but has bad type (not Int32)";
341 info() <<
"[readMeshFromVtuFile] Fetched";
347 info() <<
"[readMeshFromVtuFile] nbOfCells=" << nbOfCells <<
", nbOfNodes=" << nbOfNodes;
349 HashTableMapT<Int64,Real3> nodes_coords(nbOfNodes,
true);
350 vtkPoints* vtkAllPoints=unstructuredGrid->GetPoints();
351 for(Integer i=0; i<nbOfNodes; ++i ){
353 vtkAllPoints->GetPoint(i,xyz);
355 coords[i] = Real3(xyz[0],xyz[1],xyz[2]);
356 nodes_coords.nocheckAdd(nodesUidArray->GetValue(i),coords[i]);
360 HashTableMapT<Int64,Int32> nodes_owner_map(nbOfNodes,
true);
361 if (nodes_owner_array){
362 for(Integer i=0; i<nbOfNodes; ++i ){
363 nodes_owner_map.nocheckAdd(nodesUidArray->GetValue(i),nodes_owner_array->GetValue(i));
368 cells_filter.
resize(nbOfCells);
369 for( Integer i=0; i<nbOfCells; ++i )
374 for( Integer j=0, js=cells_filter.size(); j<js; ++j )
375 mesh_nb_cell_node += unstructuredGrid->GetCell(j)->GetNumberOfPoints();
376 info() <<
"Number of mesh_nb_cell_node = "<<mesh_nb_cell_node;
386 for( Integer i_cell=0, s_cell=cells_filter.size(); i_cell<s_cell; ++i_cell ){
387 vtkIdType iVtkCell=i_cell;
388 vtkCell *cell = unstructuredGrid->GetCell(iVtkCell);
389 if (!cell->IsLinear())
throw NotImplementedException(A_FUNCINFO);
390 Integer arcItemType = IT_NullType;
391 switch(cell->GetCellType()){
393 case(VTK_TETRA): arcItemType = IT_Tetraedron4;
break;
394 case(VTK_HEXAHEDRON): arcItemType = IT_Hexaedron8;
break;
395 case(VTK_PYRAMID): arcItemType = IT_Pyramid5;
break;
397 case(VTK_WEDGE): arcItemType = IT_Pentaedron6;
break;
398 case(VTK_PENTAGONAL_PRISM): arcItemType = IT_Heptaedron10;
break;
399 case(VTK_HEXAGONAL_PRISM): arcItemType = IT_Octaedron12;
break;
401 case(VTK_QUAD): arcItemType = IT_Quad4;
break;
402 case(VTK_TRIANGLE): arcItemType = IT_Triangle3;
break;
404 case(VTK_POLY_VERTEX):
405 switch (cell->GetNumberOfPoints()){
406 case (4): arcItemType = IT_Tetraedron4;
break;
407 case (7): arcItemType = IT_HemiHexa7;
break;
408 case (6): arcItemType = IT_HemiHexa6;
break;
409 case (5): arcItemType = IT_HemiHexa5;
break;
410 default:
throw NotImplementedException(A_FUNCINFO);;
414 switch (cell->GetNumberOfPoints()){
415 case (4): arcItemType = IT_Tetraedron4;
break;
416 default:
throw NotImplementedException(A_FUNCINFO);;
419#ifndef NO_USER_WARNING
420#warning IT_Vertex returns 0 nodes in ItemTypeMng.cc@42:type->setInfos(this,IT_Vertex,0,0,0);
421#warning IT_Vertex vs IT_DualNode HACK
423 case(VTK_VERTEX): arcItemType=IT_DualNode;
break;
425 switch (cell->GetNumberOfPoints()){
426 case (3): arcItemType = IT_Triangle3;
break;
427 case (4): arcItemType = IT_Quad4;
break;
428 case (5): arcItemType = IT_Pentagon5;
break;
429 case (6): arcItemType = IT_Hexagon6;
break;
430 default:
throw NotImplementedException(A_FUNCINFO);;
435 switch (cell->GetNumberOfPoints()){
436 case (2): arcItemType = IT_CellLine2;
break;
437 default:
throw NotImplementedException(A_FUNCINFO);;
440 case(VTK_EMPTY_CELL):
441 info() <<
"VTK_EMPTY_CELL n="<<cell->GetNumberOfPoints();
443 case(VTK_TRIANGLE_STRIP):
444 info() <<
"VTK_TRIANGLE_STRIP n="<<cell->GetNumberOfPoints();
447 info() <<
"VTK_VOXEL n="<<cell->GetNumberOfPoints();
448 switch (cell->GetNumberOfPoints()){
449 case (3): arcItemType = IT_Triangle3;
break;
450 case (4): arcItemType = IT_Quad4;
break;
451 case (5): arcItemType = IT_Pentagon5;
break;
452 case (6): arcItemType = IT_Hexagon6;
break;
453 case (7): arcItemType = IT_HemiHexa7;
break;
454 case (8): arcItemType = IT_Hexaedron8;
break;
455 default:
throw NotImplementedException(A_FUNCINFO);;
458 default:
throw NotImplementedException(A_FUNCINFO);;
460 int nNodes = cell->GetNumberOfPoints();
463 cells_infos[cells_infos_index++]=arcItemType;
466 Integer cell_indirect_id = cells_filter[i_cell];
467 cells_infos[cells_infos_index++] = cellsUidArray->GetValue(cell_indirect_id);
470 vtkIdList* nodeIds=cell->GetPointIds();
471 for(
int iNode=0; iNode<nNodes; ++iNode){
472 Integer localId = nodeIds->GetId(iNode);
473 long uniqueUid=nodesUidArray->GetValue(localId);
475 cells_infos[cells_infos_index++] = uniqueUid;
483 info() <<
"[readMeshFromVtuFile] ## Mesh 3D ##";
484 PRIMARYMESH_CAST(mesh)->setDimension(3);
485 info() <<
"[readMeshFromVtuFile] ## Allocating " << cells_filter.size() <<
" cells ##";
486 PRIMARYMESH_CAST(mesh)->allocateCells(cells_filter.size(), cells_infos,
false);
491 info() <<
"[readMeshFromVtuFile] internalNodes.size()="<<internalNodes.size();
492 if (nodes_owner_array){
493 info() <<
"Set nodes owners from vtu file";
494 for(Integer i=0, is=internalNodes.size(); i<is; ++i){
495 ItemInternal* internal_node = internalNodes[i];
496 Int32 true_owner = nodes_owner_map[internal_node->uniqueId()];
497 internal_node->setOwner(true_owner,sid);
501 for(Integer i=0, is=internalNodes.size(); i<is; ++i)
502 internalNodes[i]->setOwner(sid,sid);
505 info() <<
"[readMeshFromVtuFile] internalCells.size()="<<internalCells.size();
506 for(Integer i=0, is=internalCells.size(); i<is; ++i)
507 internalCells[i]->setOwner(sid,sid);
513 info() <<
"[readMeshFromVtuFile] ## Ending with endAllocate ##";
514 PRIMARYMESH_CAST(mesh)->endAllocate();
517 mesh->utilities()->changeOwnersFromCells();
520 info() <<
"\n[readMeshFromVtuFile] ## Now dealing with ghost's layer ##";
521 info() <<
"[readMeshFromVtuFile] mesh.nbNode=" <<mesh->nbNode() <<
" mesh.nbCell="<< mesh->nbCell();
527 info() <<
"[readMeshFromVtuFile] ## Now Fetching Other Fields ##";
528 vtkFieldData* allFieldData=unstructuredGrid->GetFieldData();
529 int nbrOfArrays = allFieldData->GetNumberOfArrays();
530 info() <<
"[readMeshFromVtuFile] nbrOfArrays = " << nbrOfArrays;
531 for(
int i=0;i<nbrOfArrays;++i){
532 if (itWasAnArcanProduction==
false)
continue;
533 info() <<
"[readMeshFromVtuFile] Focussing on \"" << allFieldData->GetArrayName(i) <<
"\" (i="<<i<<
")";
534 if (readGroupsFromFieldData(mesh, allFieldData, i) !=
true)
542 info() <<
"[readMeshFromVtuFile] ## Now insert coords ##";
546 nodes_coord_var[inode] = nodes_coords[inode->uniqueId()];
553 mesh->synchronizeGroupsAndVariables();
559 info() <<
"[readMeshFromVtuFile] RTOk";
563 nodesUidArray->Delete();
565 cellsUidArray->Delete();
#define ARCANE_REGISTER_SUB_DOMAIN_FACTORY(aclass, ainterface, aname)
Enregistre un service de fabrique pour la classe aclass.
Classe de base d'un service.
Interface du service gérant la lecture d'un maillage.
eReturnType
Types des codes de retour d'une lecture ou écriture.
@ RTError
Erreur lors de l'opération.
@ RTOk
Opération effectuée avec succès.
virtual Int32 commRank() const =0
Rang de cette instance dans le communicateur.
virtual bool isParallel() const =0
Retourne true si l'exécution est parallèle.
Interface du gestionnaire d'un sous-domaine.
Lecteur des fichiers de maillage via la bibliothèque LIMA.
Structure contenant les informations pour créer un service.
Classe utilitaire pour instantier un service d'une interface donnée.
Lecteur des fichiers de maillage aux format Vtk.
virtual eReturnType readMeshFromFile(IPrimaryMesh *mesh, const XmlNode &mesh_node, const String &file_name, const String &dir_name, bool use_internal_partition)
Lit un maillage à partir d'un fichier.
bool allowExtension(const String &str)
Vérifie si le service supporte les fichiers avec l'extension str.
virtual void build()
Construction de niveau build du service.
void resize(Int64 s)
Change le nombre d'éléments du tableau à s.
Chaîne de caractères unicode.
const char * localstr() const
Retourne la conversion de l'instance dans l'encodage UTF-8.
TraceMessage pfatal() const
Flot pour un message d'erreur fatale en parallèle.
TraceMessage warning() const
Flot pour un message d'avertissement.
TraceMessage info() const
Flot pour un message d'information.
MeshVariableScalarRefT< Node, Real3 > VariableNodeReal3
Grandeur au noeud de type coordonnées.
UniqueArray< Int64 > Int64UniqueArray
Tableau dynamique à une dimension d'entiers 64 bits.
UniqueArray< Real3 > Real3UniqueArray
Tableau dynamique à une dimension de vecteurs de rang 3.
ConstArrayView< ItemInternal * > ItemInternalList
Type de la liste interne des entités.
UniqueArray< Int32 > Int32UniqueArray
Tableau dynamique à une dimension d'entiers 32 bits.
eItemKind
Genre d'entité de maillage.
@ IK_Node
Entité de maillage de genre noeud.
@ IK_Cell
Entité de maillage de genre maille.
UniqueArray< Integer > IntegerUniqueArray
Tableau dynamique à une dimension d'entiers.
Int32 Integer
Type représentant un entier.