Arcane  v3.14.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
VtuMeshReader.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2022 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-2013 */
9/* */
10/* Lecture/Ecriture d'un fichier au format VtuMeshReader. */
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/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"
50
51#include <vtkXMLUnstructuredGridReader.h>
52#include <vtkUnstructuredGrid.h>
53#include <vtkCell.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>
60
61
62/**********************************************************************
63 * Convention data types/file type and particular file extension parity *
64 **********************************************************************/
65#define VTK_FILE_EXT_VTI "vti" // Serial structured vtkImageData
66#define VTK_FILE_EXT_VTP "vtp" // Serial UNstructured vtkPolyData
67#define VTK_FILE_EXT_VTR "vtr" // Serial structured vtkRectilinearGrid
68#define VTK_FILE_EXT_VTS "vts" // Serial structured vtkStructuredGrid
69#define VTK_FILE_EXT_VTU "vtu" // Serial UNstructured vtkUnstructuredGrid
70#define VTK_FILE_EXT_PVTI "pvti" // Parallel structured vtkImageData
71#define VTK_FILE_EXT_PVTP "pvtp" // Parallel UNstructured vtkPolyData
72#define VTK_FILE_EXT_PVTR "pvtr" // Parallel structured vtkRectilinearGrid
73#define VTK_FILE_EXT_PVTS "pvts" // Parallel structured vtkStructuredGrid
74#define VTK_FILE_EXT_PVTU "pvtu" // Parallel UNstructured vtkUnstructuredGrid
75
76
77/*---------------------------------------------------------------------------*/
78/*---------------------------------------------------------------------------*/
79
80ARCANE_BEGIN_NAMESPACE
81
82/*---------------------------------------------------------------------------*/
83/*---------------------------------------------------------------------------*/
84
85
86/*---------------------------------------------------------------------------*/
87/*---------------------------------------------------------------------------*/
92: public AbstractService, public IMeshReader
93{
94 public:
95
97
98 public:
99
100 virtual void build() {}
101
102 bool allowExtension(const String& str);
103
104 virtual eReturnType readMeshFromFile(IPrimaryMesh* mesh,const XmlNode& mesh_node,
106 eReturnType readMeshFromVtuFile(IMesh* mesh,const XmlNode& mesh_node,
108
109 ISubDomain* subDomain() { return m_sub_domain; }
110
111 bool readGroupsFromFieldData(IMesh *mesh, vtkFieldData*, int);
112
113 private:
114
115 int vtkFileExtIdx;
116 ISubDomain* m_sub_domain;
117 eReturnType writeMeshToLemFile(IMesh* mesh, const String& filename,const String& dir_name);
118};
119
120/*---------------------------------------------------------------------------*/
121/*---------------------------------------------------------------------------*/
122
123ARCANE_REGISTER_SUB_DOMAIN_FACTORY(VtuMeshReader,IMeshReader,VtuNewMeshReader);
124
125/*---------------------------------------------------------------------------*/
126/*---------------------------------------------------------------------------*/
127
128/****************************************************************************
129 * Assignation des data types/file-type and particular file extension readers *
130\****************************************************************************/
131typedef struct{
132 char *ext;
133 IMeshReader::eReturnType (VtuMeshReader::*reader)(IMesh*,const XmlNode&,const String&,const String&,bool);
135
136vtkExtReader vtkFileExtReader[]={
137 {VTK_FILE_EXT_VTU, &VtuMeshReader::readMeshFromVtuFile},
138 {NULL,NULL}
139};
140
141
142/*---------------------------------------------------------------------------*/
143/*---------------------------------------------------------------------------*/
144
145VtuMeshReader::
146VtuMeshReader(const ServiceBuildInfo& sbi)
147: AbstractService(sbi)
148, m_sub_domain(sbi.subDomain())
149{
150 vtkFileExtIdx=0;
151}
152
153/*---------------------------------------------------------------------------*/
154/*---------------------------------------------------------------------------*/
155
156/*****************************************************************
157 * allowExtension
158 *****************************************************************/
160allowExtension(const String& str)
161{
162 //info() << "[allowExtension] Checking for file extension...";
163 for(vtkFileExtIdx=0;vtkFileExtReader[vtkFileExtIdx].ext!=NULL;++vtkFileExtIdx){
164 //info() << "Testing for '" << vtkFileExtReader[vtkFileExtIdx].ext << "'...";
165 if (str == vtkFileExtReader[vtkFileExtIdx].ext){
166 return true;
167 }
168 }
169 //info()<<"Miss for all!";
170 // Sets our index in place for further service management
171 vtkFileExtIdx=0;
172 return false;
173}
174
175/*---------------------------------------------------------------------------*/
176/*---------------------------------------------------------------------------*/
177
178/*****************************************************************
179 * readMeshFromFile switch
180 *****************************************************************/
183 const String& filename,const String& dir_name,
185{
186 info() << "[readMeshFromFile] Forwarding to vtkFileExtReader[" << vtkFileExtIdx << "].reader";
187 return (this->*vtkFileExtReader[vtkFileExtIdx].reader)(mesh, mesh_node, filename, dir_name, use_internal_partition);
188}
189
190/*---------------------------------------------------------------------------*/
191/*---------------------------------------------------------------------------*/
192
193/*****************************************************************
194 * writeMeshToLemFile
195 *****************************************************************/
196IMeshReader::eReturnType VtuMeshReader::
197writeMeshToLemFile(IMesh* mesh, const String& filename,const String& dir_name)
198{
199 ISubDomain* sd = m_sub_domain;
201 //FactoryT<IMeshWriter> mesh_writer_factory(sd->serviceMng());
202 //mesh_writer = mesh_writer_factory.createInstance("Lima",true);
203 if (!mesh_writer.get())
204 pfatal() << "Mesh writer service selected not available";
205 std::string fname = filename.localstr();
206 fname += ".unf";
207 if (mesh_writer->writeMeshToFile(mesh,fname) != RTOk)
208 return RTError;
209 return RTOk;
210}
211
212/*---------------------------------------------------------------------------*/
213/*---------------------------------------------------------------------------*/
214/*****************************************************************
215 * readGroupsFromFieldData
216 *****************************************************************/
217bool VtuMeshReader::
218readGroupsFromFieldData(IMesh* mesh,vtkFieldData* allFieldData,int i)
219{
220 const std::string group_name = allFieldData->GetArrayName(i);
221 vtkDataArray* iFieldData = allFieldData->GetArray(i);
222 if (!iFieldData)
223 return false;
224
225 if (iFieldData->GetDataType() != VTK_LONG)
226 return false;
227
228 Integer nb_tuple = iFieldData->GetNumberOfTuples();
229 info() << "[readGroupsFromFieldData] iFieldData->GetNumberOfTuples="<< nb_tuple;
230 if (nb_tuple==0)
231 return false;
232 vtkLongArray* vtk_array = vtkLongArray::SafeDownCast(iFieldData);
233
234 // Le premier élément du tableau contient son type
235 eItemKind kind_type = (eItemKind)(vtk_array->GetValue(0));
236 IItemFamily* family = mesh->itemFamily(kind_type);
237 Integer nb_item = nb_tuple - 1;
238 Int64UniqueArray unique_ids(nb_item);
239 // Les éléments suivant contiennent les uniqueId() des entités du groupe.
240 for( Integer z=0; z<nb_item; ++z )
241 unique_ids[z] = vtk_array->GetValue(z+1);
242
243 // Récupère le localId() correspondant.
244 Int32UniqueArray local_ids(unique_ids.size());
245 family->itemsUniqueIdToLocalId(local_ids,unique_ids,false);
246
247 // Tous les entités ne sont pas forcément dans le maillage actuel et
248 // il faut donc les filtrer.
250 for( Integer i=0; i<nb_item; ++i )
251 if (local_ids[i]!=NULL_ITEM_LOCAL_ID)
252 ids.add(local_ids[i]);
253
254 info() << "Create group family=" << family->name() << " name=" << group_name << " ids=" << ids.size();
255 family->createGroup(group_name,ids);
256 return true;
257}
258
259/*---------------------------------------------------------------------------*/
260/*---------------------------------------------------------------------------*/
261
262/*****************************************************************
263 * readMeshFromVtuFile
264 *****************************************************************/
265IMeshReader::eReturnType VtuMeshReader::
266readMeshFromVtuFile(IMesh* mesh, const XmlNode& mesh_node,
267 const String& filename,const String& dir_name,
268 bool use_internal_partition)
269{
270 bool itWasAnArcanProduction=true;
271 IParallelMng* pm = mesh->parallelMng();
272 bool is_parallel = pm->isParallel();
273 Integer sid = pm->commRank();
274
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()))
280 return RTError;
281 reader->UpdateInformation();
282 reader->Update();// Force reading
283 vtkUnstructuredGrid *unstructuredGrid = reader->GetOutput();
284 // Avec VTK 7, plus besoin du Update.
285 //unstructuredGrid->Update();// La lecture effective du fichier n'a lieu qu'après l'appel à Update().
286 int nbOfCells = unstructuredGrid->GetNumberOfCells();
287 int nbOfNodes = unstructuredGrid->GetNumberOfPoints();
288
289
290 /*******************************
291 *Fetching Nodes UID & Cells UID *
292 *******************************/
293 info() << "[readMeshFromVtuFile] ## Now Fetching Nodes Unique IDs ##";
294 vtkPointData* allPointData=unstructuredGrid->GetPointData(); // Data associated to Points
295 allPointData->Update();
296 vtkDataArray* dataNodeArray = allPointData->GetArray("NodesUniqueIDs");
297 vtkLongArray *nodesUidArray = nullptr;
298 if (!dataNodeArray){
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;
304 }
305 else{
306 if (dataNodeArray->GetDataType() != VTK_LONG)
307 return RTError;
308 nodesUidArray = vtkLongArray::SafeDownCast(dataNodeArray);
309 }
310 info() << "[readMeshFromVtuFile] Fetched";
311
312 info() << "[readMeshFromVtuFile] ## Now Fetching Cells Unique IDs ##";
313 vtkCellData* allCellData=unstructuredGrid->GetCellData(); // Data associated to Points
314 allCellData->Update();
315 vtkDataArray* dataCellArray = allCellData->GetArray("CellsUniqueIDs");
316 vtkLongArray *cellsUidArray = nullptr;
317 if (!dataCellArray){
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;
323 }
324 else{
325 if (dataCellArray->GetDataType() != VTK_LONG)
326 return RTError;
327 cellsUidArray = vtkLongArray::SafeDownCast(dataCellArray);
328 }
329
330 // Tableau contenant les numéros des propriétaires des noeuds.
331 // Ce tableau est optionnel et est utilisé par le partitionneur
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);
337 }
338 else
339 warning() << "NodesOwner array is present but has bad type (not Int32)";
340 }
341 info() << "[readMeshFromVtuFile] Fetched";
342
343
344 /************************
345 * Fetch own nodes coords *
346 ************************/
347 info() << "[readMeshFromVtuFile] nbOfCells=" << nbOfCells << ", nbOfNodes=" << nbOfNodes;
348 Real3UniqueArray coords(nbOfNodes);
349 HashTableMapT<Int64,Real3> nodes_coords(nbOfNodes,true);
350 vtkPoints* vtkAllPoints=unstructuredGrid->GetPoints();
351 for(Integer i=0; i<nbOfNodes; ++i ){
352 double xyz[3];
353 vtkAllPoints->GetPoint(i,xyz);
354 //info() << "x=" << xyz[0] << " y=" << xyz[1]<< " z=" << xyz[2] << ", nodUid=" << nodesUidArray->GetValue(i);
355 coords[i] = Real3(xyz[0],xyz[1],xyz[2]);
356 nodes_coords.nocheckAdd(nodesUidArray->GetValue(i),coords[i]);
357 }
358
359 // Create hash table for nodes owner.
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));
364 }
365 }
366
367 IntegerUniqueArray cells_filter;
368 cells_filter.resize(nbOfCells);
369 for( Integer i=0; i<nbOfCells; ++i )
370 cells_filter[i] = i;
371
372 // Calcul le nombre de mailles/noeuds
373 Integer mesh_nb_cell_node = 0;
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;
377
378 // Tableau contenant les infos aux mailles (voir IMesh::allocateMesh())
379 Int64UniqueArray cells_infos(mesh_nb_cell_node+cells_filter.size()*2);
380 Integer cells_infos_index = 0;
381
382 /******************
383 * Now insert CELLS *
384 ******************/
385 // For all cells that are discovered
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()){
392 // Linear cells
393 case(VTK_TETRA): arcItemType = IT_Tetraedron4; break;
394 case(VTK_HEXAHEDRON): arcItemType = IT_Hexaedron8; break;
395 case(VTK_PYRAMID): arcItemType = IT_Pyramid5; break;
396 /* Others not yet implemented */
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;
400 /* 2D */
401 case(VTK_QUAD): arcItemType = IT_Quad4; break;
402 case(VTK_TRIANGLE): arcItemType = IT_Triangle3; break;
403 // Cas du poly vertex à parser
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);;
411 }
412 break;
413 case(VTK_PIXEL):
414 switch (cell->GetNumberOfPoints()){
415 case (4): arcItemType = IT_Tetraedron4; break;
416 default:throw NotImplementedException(A_FUNCINFO);;
417 }
418 break;
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
422#endif /* NO_USER_WARNING */
423 case(VTK_VERTEX): arcItemType=IT_DualNode; break;
424 case(VTK_POLYGON):
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);;
431 }
432 break;
433 case(VTK_LINE):
434 case(VTK_POLY_LINE):
435 switch (cell->GetNumberOfPoints()){
436 case (2): arcItemType = IT_CellLine2; break;
437 default:throw NotImplementedException(A_FUNCINFO);;
438 }
439 break;
440 case(VTK_EMPTY_CELL):
441 info() << "VTK_EMPTY_CELL n="<<cell->GetNumberOfPoints();
442 break;
443 case(VTK_TRIANGLE_STRIP):
444 info() << "VTK_TRIANGLE_STRIP n="<<cell->GetNumberOfPoints();
445 break;
446 case(VTK_VOXEL):
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);;
456 }
457 break;
458 default:throw NotImplementedException(A_FUNCINFO);;
459 }
460 int nNodes = cell->GetNumberOfPoints();// Return the number of points in the cell
461
462 // First is cell's TYPE Stocke le type de la maille
463 cells_infos[cells_infos_index++]=arcItemType;
464
465 // Then comes its UniqueID Stocke le numéro unique de la maille
466 Integer cell_indirect_id = cells_filter[i_cell];
467 cells_infos[cells_infos_index++] = cellsUidArray->GetValue(cell_indirect_id);
468
469 // And finally the Nodes' unique IDs
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);
474 //info() << "working on localId=" << localId << ", uniqueUid=" << uniqueUid;
475 cells_infos[cells_infos_index++] = uniqueUid;
476 }
477 }
478
479
480 /********************************
481 * Setting Dimension & Allocating *
482 ********************************/
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);
487
488
489 // Positionne les propriétaires des noeuds à partir des groupes de noeuds
490 ItemInternalList internalNodes(mesh->itemsInternal(IK_Node));
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);
498 }
499 }
500 else
501 for(Integer i=0, is=internalNodes.size(); i<is; ++i)
502 internalNodes[i]->setOwner(sid,sid);
503
504 ItemInternalList internalCells(mesh->itemsInternal(IK_Cell));
505 info() << "[readMeshFromVtuFile] internalCells.size()="<<internalCells.size();
506 for(Integer i=0, is=internalCells.size(); i<is; ++i)
507 internalCells[i]->setOwner(sid,sid);
508
509
510 /********************************************
511 * Now finishing & preparing for ghost layout *
512 ********************************************/
513 info() << "[readMeshFromVtuFile] ## Ending with endAllocate ##";
514 PRIMARYMESH_CAST(mesh)->endAllocate();
515 if (is_parallel) {
516 // mesh->setOwnersFromCells();
517 mesh->utilities()->changeOwnersFromCells();
518 }
519
520 info() << "\n[readMeshFromVtuFile] ## Now dealing with ghost's layer ##";
521 info() << "[readMeshFromVtuFile] mesh.nbNode=" <<mesh->nbNode() << " mesh.nbCell="<< mesh->nbCell();
522
523
524 /***********************
525 * Fetching Other Groups *
526 ***********************/
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)
535 return RTError;
536 }
537
538
539 /*******************
540 * Now insert coords *
541 *******************/
542 info() << "[readMeshFromVtuFile] ## Now insert coords ##";
543 // Remplit la variable contenant les coordonnées des noeuds
544 VariableNodeReal3& nodes_coord_var(PRIMARYMESH_CAST(mesh)->nodesCoordinates());
545 ENUMERATE_NODE(inode,mesh->ownNodes()){
546 nodes_coord_var[inode] = nodes_coords[inode->uniqueId()];
547 }
548
549
550 /****************************************
551 * Synchronizing groups/variables & nodes *
552 ****************************************/
553 mesh->synchronizeGroupsAndVariables();
554
555 /**************
556 * Finishing up *
557 **************/
558 reader->Delete();
559 info() << "[readMeshFromVtuFile] RTOk";
560
561 // TODO: regarder comment detruire automatiquement
562 if (!dataNodeArray)
563 nodesUidArray->Delete();
564 if (!dataCellArray)
565 cellsUidArray->Delete();
566
567 return RTOk;
568}
569
570/*---------------------------------------------------------------------------*/
571/*---------------------------------------------------------------------------*/
572
573ARCANE_END_NAMESPACE
574
575/*---------------------------------------------------------------------------*/
576/*---------------------------------------------------------------------------*/
#define ENUMERATE_NODE(name, group)
Enumérateur générique d'un groupe de noeuds.
#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.
Definition IMeshReader.h:37
eReturnType
Types des codes de retour d'une lecture ou écriture.
Definition IMeshReader.h:42
@ RTError
Erreur lors de l'opération.
Definition IMeshReader.h:44
@ RTOk
Opération effectuée avec succès.
Definition IMeshReader.h:43
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.
Definition ISubDomain.h:74
Lecteur des fichiers de maillage via la bibliothèque LIMA.
Definition Lima.cc:120
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.
Noeud d'un arbre DOM.
Definition XmlNode.h:51
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.
Definition String.cc:227
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.
Definition UtilsTypes.h:513
UniqueArray< Real3 > Real3UniqueArray
Tableau dynamique à une dimension de vecteurs de rang 3.
Definition UtilsTypes.h:529
ConstArrayView< ItemInternal * > ItemInternalList
Type de la liste interne des entités.
Definition ItemTypes.h:466
UniqueArray< Int32 > Int32UniqueArray
Tableau dynamique à une dimension d'entiers 32 bits.
Definition UtilsTypes.h:515
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.
Definition UtilsTypes.h:519
Int32 Integer
Type représentant un entier.