Arcane  v3.15.0.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
VtkMeshIOService.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2024 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/* VtkMeshIOService.cc (C) 2000-2024 */
9/* */
10/* Lecture/Ecriture d'un maillage au format Vtk historique (legacy). */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
14#include "arcane/utils/Collection.h"
15#include "arcane/utils/Enumerator.h"
16#include "arcane/utils/HashTableMap.h"
17#include "arcane/utils/IOException.h"
18#include "arcane/utils/ITraceMng.h"
19#include "arcane/utils/Iostream.h"
20#include "arcane/utils/OStringStream.h"
21#include "arcane/utils/ScopedPtr.h"
22#include "arcane/utils/StdHeader.h"
23#include "arcane/utils/String.h"
24#include "arcane/utils/ValueConvert.h"
25#include "arcane/utils/Real3.h"
26
27#include "arcane/core/BasicService.h"
28#include "arcane/core/FactoryService.h"
29#include "arcane/core/ICaseMeshReader.h"
30#include "arcane/core/IItemFamily.h"
31#include "arcane/core/IPrimaryMesh.h"
32#include "arcane/core/IMeshBuilder.h"
33#include "arcane/core/IMeshReader.h"
34#include "arcane/core/IMeshUtilities.h"
35#include "arcane/core/IMeshWriter.h"
36#include "arcane/core/IParallelMng.h"
37#include "arcane/core/IVariableAccessor.h"
38#include "arcane/core/IXmlDocumentHolder.h"
39#include "arcane/core/Item.h"
41#include "arcane/core/IVariableMng.h"
42#include "arcane/core/VariableTypes.h"
43#include "arcane/core/XmlNode.h"
44#include "arcane/core/XmlNodeList.h"
45#include "arcane/core/UnstructuredMeshAllocateBuildInfo.h"
46
47#include "arcane/core/internal/IVariableMngInternal.h"
48#include "arcane/core/internal/VtkCellTypes.h"
49
50/*---------------------------------------------------------------------------*/
51/*---------------------------------------------------------------------------*/
52
53namespace Arcane
54{
55
56/*---------------------------------------------------------------------------*/
57/*---------------------------------------------------------------------------*/
58
59class VtkFile;
60using namespace Arcane::VtkUtils;
61
62/*---------------------------------------------------------------------------*/
63/*---------------------------------------------------------------------------*/
102: public TraceAccessor
103{
104 public:
105
106 explicit VtkMeshIOService(ITraceMng* tm)
108 {}
109
110 public:
111
112 void build() {}
113
114 public:
115
116 enum eMeshType
117 {
118 VTK_MT_Unknown,
119 VTK_MT_StructuredGrid,
120 VTK_MT_UnstructuredGrid
121 };
122
124 {
125 public:
126 };
127
129 : public VtkMesh
130 {
131 public:
132
133 int m_nb_x;
134 int m_nb_y;
135 int m_nb_z;
136 };
137
138 public:
139
141
142 private:
143
146 void _readCellVariable(IMesh* mesh, VtkFile& vtk_file, const String& name_str, Integer nb_cell);
147 void _readItemGroup(IMesh* mesh, VtkFile& vtk_file, const String& name_str, Integer nb_item,
149 void _readNodeGroup(IMesh* mesh, VtkFile& vtk_file, const String& name, Integer nb_item);
152 Int32ConstArrayView local_id, Integer nb_node);
155 Array<Integer>& cells_nb_node,
156 Array<ItemTypeId>& cells_type,
157 Array<Int64>& cells_connectivity);
158 void _readFacesMesh(IMesh* mesh, const String& file_name,
160 bool _readMetadata(IMesh* mesh, VtkFile& vtk_file);
161
162 private:
163};
164
165/*---------------------------------------------------------------------------*/
166/*---------------------------------------------------------------------------*/
167
169{
170 public:
171
172 static const int BUFSIZE = 10000;
173
174 public:
175
176 explicit VtkFile(std::istream* stream)
177 : m_stream(stream)
178 , m_is_init(false)
180 , m_is_eof(false)
181 , m_is_binary_file(false)
182 , m_buf{}
183 {}
184
185 const char* getCurrentLine();
186 bool isEmptyNextLine();
187 const char* getNextLine();
188
190 void checkString(const String& current_value,
191 const String& expected_value1,
192 const String& expected_value2);
193
194 static bool isEqualString(const String& current_value, const String& expected_value);
195
196 void reReadSameLine() { m_need_reread_current_line = true; }
197
198 bool isEof() { return m_is_eof; }
199
200 template <class T>
201 void getBinary(T& type);
202 float getFloat();
203 double getDouble();
204 int getInt();
205
206 void setIsBinaryFile(bool new_val) { m_is_binary_file = new_val; }
207
208 private:
209
211 std::istream* m_stream = nullptr;
212
215
218
221
224
226 char m_buf[BUFSIZE];
227};
228
229/*---------------------------------------------------------------------------*/
230/*---------------------------------------------------------------------------*/
231
237const char* VtkFile::
239{
240 if (!m_is_init)
241 getNextLine();
242 return m_buf;
243}
244
245/*---------------------------------------------------------------------------*/
246/*---------------------------------------------------------------------------*/
247
259{
260 m_is_init = true;
261
262 // On veut que getNextLine lise une nouvelle ligne.
263 // (on met à false dans le cas où cette méthode serai
264 // appelée plusieurs fois à la suite).
266
267 // Si l'on est arrivé à la fin du fichier lors du précédent appel de cette méthode ou
268 // de getNextLine, on envoie une erreur.
269 if (m_is_eof) {
270 throw IOException("VtkFile::isEmptyNextLine()", "Unexpected EndOfFile");
271 }
272
273 if (m_stream->good()) {
274 // Le getline s'arrete (par défaut) au char '\n' et ne l'inclus pas dans le buf
275 // mais le remplace par '\0'.
276 m_stream->getline(m_buf, sizeof(m_buf) - 1);
277
278 // Si on arrive au bout du fichier, on return true (pour dire oui, il y a une ligne vide,
279 // à l'appelant de gérer ça).
280 if (m_stream->eof()) {
281 m_is_eof = true;
282 return true;
283 }
284
285 // Sous Windows, une ligne vide commence par \r.
286 // getline remplace \n par \0, que ce soit sous Windows ou Linux.
287 if (m_buf[0] == '\r' || m_buf[0] == '\0') {
288 getNextLine();
289
290 // On demande à ce que le prochain appel à getNextLine renvoie la ligne
291 // qui vient tout juste d'être bufferisée.
293 return true;
294 }
295 else {
296 bool is_comment = true;
297
298 // On retire le commentaire, s'il y en a un, en remplaçant '#' par '\0'.
299 for (int i = 0; i < BUFSIZE && m_buf[i] != '\0'; ++i) {
300 if (!isspace(m_buf[i]) && m_buf[i] != '#' && is_comment) {
301 is_comment = false;
302 }
303 if (m_buf[i] == '#') {
304 m_buf[i] = '\0';
305 break;
306 }
307 }
308
309 // Si ce n'est pas un commentaire, on supprime juste le '\r' final (si windows).
310 if (!is_comment) {
311 // Supprime le '\r' final
312 for (int i = 0; i < BUFSIZE && m_buf[i] != '\0'; ++i) {
313 if (m_buf[i] == '\r') {
314 m_buf[i] = '\0';
315 break;
316 }
317 }
318 }
319
320 // Si c'était un commentaire, on recherche la prochaine ligne "valide"
321 // en appelant getNextLine.
322 else {
323 getNextLine();
324 }
325 }
327 return false;
328 }
329 throw IOException("VtkFile::isEmptyNextLine()", "Not Good");
330}
331
332/*---------------------------------------------------------------------------*/
333/*---------------------------------------------------------------------------*/
334
340const char* VtkFile::
342{
343 m_is_init = true;
344
345 // On return le buffer actuel, si celui-ci n'a pas été utilisé.
348 return getCurrentLine();
349 }
350
351 // Si l'on est arrivé à la fin du fichier lors du précédent appel de cette méthode ou
352 // de isEmptyNextLine, on envoie une erreur.
353 if (m_is_eof) {
354 throw IOException("VtkFile::isEmptyNextLine()", "Unexpected EndOfFile");
355 }
356
357 while (m_stream->good()) {
358 // Le getline s'arrete (par défaut) au char '\n' et ne l'inclus pas dans le buf mais le remplace par '\0'.
359 m_stream->getline(m_buf, sizeof(m_buf) - 1);
360
361 // Si on arrive au bout du fichier, on return le buffer avec \0 au début (c'est à l'appelant d'appeler
362 // isEof() pour savoir si le fichier est fini ou non).
363 if (m_stream->eof()) {
364 m_is_eof = true;
365 m_buf[0] = '\0';
366 return m_buf;
367 }
368
369 bool is_comment = true;
370
371 // Sous Windows, une ligne vide commence par \r.
372 // getline remplace \n par \0, que ce soit sous Windows ou Linux.
373 if (m_buf[0] == '\0' || m_buf[0] == '\r')
374 continue;
375
376 // On retire le commentaire, s'il y en a un, en remplaçant '#' par '\0'.
377 for (int i = 0; i < BUFSIZE && m_buf[i] != '\0'; ++i) {
378 if (!isspace(m_buf[i]) && m_buf[i] != '#' && is_comment) {
379 is_comment = false;
380 }
381 if (m_buf[i] == '#') {
382 m_buf[i] = '\0';
383 break;
384 }
385 }
386
387 // Si ce n'est pas un commentaire, on supprime juste le '\r' final (si windows).
388 if (!is_comment) {
389 // Supprime le '\r' final
390 for (int i = 0; i < BUFSIZE && m_buf[i] != '\0'; ++i) {
391 if (m_buf[i] == '\r') {
392 m_buf[i] = '\0';
393 break;
394 }
395 }
396 return m_buf;
397 }
398 }
399 throw IOException("VtkFile::getNextLine()", "Not good");
400}
401
402/*---------------------------------------------------------------------------*/
403/*---------------------------------------------------------------------------*/
404
411getFloat()
412{
413 float v = 0.;
414 if (m_is_binary_file) {
415 getBinary(v);
416 return v;
417 }
418 (*m_stream) >> ws >> v;
419
420 if (m_stream->good())
421 return v;
422
423 throw IOException("VtkFile::getFloat()", "Bad float");
424}
425
426/*---------------------------------------------------------------------------*/
427/*---------------------------------------------------------------------------*/
428
435getDouble()
436{
437 double v = 0.;
438 if (m_is_binary_file) {
439 getBinary(v);
440 return v;
441 }
442 (*m_stream) >> ws >> v;
443
444 if (m_stream->good())
445 return v;
446
447 throw IOException("VtkFile::getDouble()", "Bad double");
448}
449
450/*---------------------------------------------------------------------------*/
451/*---------------------------------------------------------------------------*/
452
459getInt()
460{
461 int v = 0;
462 if (m_is_binary_file) {
463 getBinary(v);
464 return v;
465 }
466 (*m_stream) >> ws >> v;
467
468 if (m_stream->good())
469 return v;
470
471 throw IOException("VtkFile::getInt()", "Bad int");
472}
473
474/*---------------------------------------------------------------------------*/
475/*---------------------------------------------------------------------------*/
476
482template <class T>
484getBinary(T& value)
485{
486 constexpr size_t sizeofT = sizeof(T);
487
488 // Le fichier VTK est en big endian et les CPU actuels sont en little endian.
491
492 // On lit les 'sizeofT' prochains octets que l'on met dans big_endian.
493 m_stream->read((char*)big_endian, sizeofT);
494
495 // On transforme le big_endian en little_endian.
496 for (size_t i = 0; i < sizeofT; i++) {
497 little_endian[sizeofT - 1 - i] = big_endian[i];
498 }
499
500 // On 'cast' la liste d'octet en type 'T'.
501 T* conv = new (little_endian) T;
502 value = *conv;
503}
504
505/*---------------------------------------------------------------------------*/
506/*---------------------------------------------------------------------------*/
507
519{
522
524 String s = "Expecting chain '" + expected_value + "', found '" + current_value + "'";
525 throw IOException("VtkFile::checkString()", s);
526 }
527}
528
529/*---------------------------------------------------------------------------*/
530/*---------------------------------------------------------------------------*/
531
544{
548
550 String s = "Expecting chain '" + expected_value1 + "' or '" + expected_value2 + "', found '" + current_value + "'";
551 throw IOException("VtkFile::checkString()", s);
552 }
553}
554
555/*---------------------------------------------------------------------------*/
556/*---------------------------------------------------------------------------*/
557
574
575/*---------------------------------------------------------------------------*/
576/*---------------------------------------------------------------------------*/
577
578/*---------------------------------------------------------------------------*/
579/*---------------------------------------------------------------------------*/
580
592{
593 std::ifstream ifile(file_name.localstr(), std::ifstream::binary);
594
595 if (!ifile) {
596 error() << "Unable to read file '" << file_name << "'";
597 return true;
598 }
599
600 debug() << "Fichier ouvert : " << file_name.localstr();
601
603 const char* buf = 0;
604
605 // Lecture de la description
606 // Lecture title.
607 String title = vtk_file.getNextLine();
608
609 info() << "Titre du fichier VTK : " << title.localstr();
610
611 // Lecture format.
612 String format = vtk_file.getNextLine();
613
614 debug() << "Format du fichier VTK : " << format.localstr();
615
616 if (VtkFile::isEqualString(format, "BINARY")) {
617 vtk_file.setIsBinaryFile(true);
618 }
619
620 eMeshType mesh_type = VTK_MT_Unknown;
621
622 // Lecture du type de maillage
623 // TODO: en parallèle, avec use_internal_partition vrai, seul le processeur 0
624 // lit les données. Dans ce cas, inutile que les autres ouvrent le fichier.
625 {
626 buf = vtk_file.getNextLine();
627
628 std::istringstream mesh_type_line(buf);
629 std::string dataset_str;
630 std::string mesh_type_str;
631
632 mesh_type_line >> ws >> dataset_str >> ws >> mesh_type_str;
633
634 vtk_file.checkString(dataset_str, "DATASET");
635
636 if (VtkFile::isEqualString(mesh_type_str, "STRUCTURED_GRID")) {
637 mesh_type = VTK_MT_StructuredGrid;
638 }
639
640 if (VtkFile::isEqualString(mesh_type_str, "UNSTRUCTURED_GRID")) {
641 mesh_type = VTK_MT_UnstructuredGrid;
642 }
643
644 if (mesh_type == VTK_MT_Unknown) {
645 error() << "Support exists only for 'STRUCTURED_GRID' and 'UNSTRUCTURED_GRID' formats (format=" << mesh_type_str << "')";
646 return true;
647 }
648 }
649 debug() << "Lecture en-tête OK";
650
651 bool ret = true;
652 switch (mesh_type) {
653 case VTK_MT_StructuredGrid:
655 break;
656
657 case VTK_MT_UnstructuredGrid:
659 debug() << "Lecture _readUnstructuredGrid OK";
660 if (!ret) {
661 // Tente de lire le fichier des faces s'il existe
663 debug() << "Lecture _readFacesMesh OK";
664 }
665 break;
666
667 case VTK_MT_Unknown:
668 break;
669 }
670 /*while ( (buf=vtk_file.getNextLine()) != 0 ){
671 info() << " STR " << buf;
672 }*/
673
674 ifile.close();
675 return ret;
676}
677
678/*---------------------------------------------------------------------------*/
679/*---------------------------------------------------------------------------*/
680
691{
692 // Lecture du nombre de points: DIMENSIONS nx ny nz
693 const char* buf = nullptr;
694 Integer nb_node_x = 0;
695 Integer nb_node_y = 0;
696 Integer nb_node_z = 0;
697 {
698 buf = vtk_file.getNextLine();
699 std::istringstream iline(buf);
700 std::string dimension_str;
701 iline >> ws >> dimension_str >> ws >> nb_node_x >> ws >> nb_node_y >> ws >> nb_node_z;
702
703 if (!iline) {
704 error() << "Syntax error while reading grid dimensions";
705 return true;
706 }
707
708 vtk_file.checkString(dimension_str, "DIMENSIONS");
709 if (nb_node_x <= 1 || nb_node_y <= 1 || nb_node_z <= 1) {
710 error() << "Invalid dimensions: x=" << nb_node_x << " y=" << nb_node_y << " z=" << nb_node_z;
711 return true;
712 }
713 }
714 info() << " Infos: " << nb_node_x << " " << nb_node_y << " " << nb_node_z;
715 Integer nb_node = nb_node_x * nb_node_y * nb_node_z;
716
717 // Lecture du nombre de points: POINTS nb float
718 std::string float_str;
719 {
720 buf = vtk_file.getNextLine();
721 std::istringstream iline(buf);
722 std::string points_str;
723 Integer nb_node_read = 0;
724 iline >> ws >> points_str >> ws >> nb_node_read >> ws >> float_str;
725 if (!iline) {
726 error() << "Syntax error while reading grid dimensions";
727 return true;
728 }
729 vtk_file.checkString(points_str, "POINTS");
730 if (nb_node_read != nb_node) {
731 error() << "Number of invalid nodes: expected=" << nb_node << " found=" << nb_node_read;
732 return true;
733 }
734 }
735
736 Int32 rank = mesh->parallelMng()->commRank();
737
738 Integer nb_cell_x = nb_node_x - 1;
739 Integer nb_cell_y = nb_node_y - 1;
740 Integer nb_cell_z = nb_node_z - 1;
741
742 if (use_internal_partition && rank != 0) {
743 nb_node_x = 0;
744 nb_node_y = 0;
745 nb_node_z = 0;
746 nb_cell_x = 0;
747 nb_cell_y = 0;
748 nb_cell_z = 0;
749 }
750
751 const Integer nb_node_yz = nb_node_y * nb_node_z;
752 const Integer nb_node_xy = nb_node_x * nb_node_y;
753
754 Integer nb_cell = nb_cell_x * nb_cell_y * nb_cell_z;
756
757 // Creation du maillage
758 {
759 UniqueArray<Integer> nodes_unique_id(nb_node);
760
761 info() << " NODE YZ = " << nb_node_yz;
762 // Création des noeuds
763 //Integer nb_node_local_id = 0;
764 {
765 Integer node_local_id = 0;
766 for (Integer x = 0; x < nb_node_x; ++x) {
767 for (Integer z = 0; z < nb_node_z; ++z) {
768 for (Integer y = 0; y < nb_node_y; ++y) {
769
770 Integer node_unique_id = y + (z)*nb_node_y + x * nb_node_y * nb_node_z;
771
772 nodes_unique_id[node_local_id] = node_unique_id;
773
775 }
776 }
777 }
778 //nb_node_local_id = node_local_id;
779 //warning() << " NB NODE LOCAL ID=" << node_local_id;
780 }
781
782 // Création des mailles
783
784 // Infos pour la création des mailles
785 // par maille: 1 pour son unique id,
786 // 1 pour son type,
787 // 8 pour chaque noeud
789
790 {
791 Integer cell_local_id = 0;
792 Integer cells_infos_index = 0;
793
794 // Normalement ne doit pas arriver car les valeurs de nb_node_x et
795 // nb_node_y sont testées lors de la lecture.
796 if (nb_node_xy == 0)
797 ARCANE_FATAL("Null value for nb_node_xy");
798
799 //Integer index = 0;
800 for (Integer z = 0; z < nb_cell_z; ++z) {
801 for (Integer y = 0; y < nb_cell_y; ++y) {
802 for (Integer x = 0; x < nb_cell_x; ++x) {
803 Integer current_cell_nb_node = 8;
804
805 //Integer cell_unique_id = y + (z)*nb_cell_y + x*nb_cell_y*nb_cell_z;
806 Int64 cell_unique_id = x + y * nb_cell_x + z * nb_cell_x * nb_cell_y;
807
808 cells_infos[cells_infos_index] = IT_Hexaedron8;
810
813
814 //Integer base_id = y + z*nb_node_y + x*nb_node_yz;
815 Integer base_id = x + y * nb_node_x + z * nb_node_xy;
816 cells_infos[cells_infos_index + 0] = nodes_unique_id[base_id];
817 cells_infos[cells_infos_index + 1] = nodes_unique_id[base_id + 1];
818 cells_infos[cells_infos_index + 2] = nodes_unique_id[base_id + nb_node_x + 1];
819 cells_infos[cells_infos_index + 3] = nodes_unique_id[base_id + nb_node_x + 0];
820 cells_infos[cells_infos_index + 4] = nodes_unique_id[base_id + nb_node_xy];
821 cells_infos[cells_infos_index + 5] = nodes_unique_id[base_id + nb_node_xy + 1];
822 cells_infos[cells_infos_index + 6] = nodes_unique_id[base_id + nb_node_xy + nb_node_x + 1];
823 cells_infos[cells_infos_index + 7] = nodes_unique_id[base_id + nb_node_xy + nb_node_x + 0];
827 }
828 }
829 }
830 }
831
832 mesh->setDimension(3);
833 mesh->allocateCells(nb_cell, cells_infos, false);
834 mesh->endAllocate();
835
836 // Positionne les coordonnées
837 {
839
840 if (vtk_file.isEqualString(float_str, "int")) {
841 for (Integer z = 0; z < nb_node_z; ++z) {
842 for (Integer y = 0; y < nb_node_y; ++y) {
843 for (Integer x = 0; x < nb_node_x; ++x) {
844 Real nx = vtk_file.getInt();
845 Real ny = vtk_file.getInt();
846 Real nz = vtk_file.getInt();
847 Integer node_unique_id = x + y * nb_node_x + z * nb_node_xy;
848 coords[node_unique_id] = Real3(nx, ny, nz);
849 }
850 }
851 }
852 }
853 else if (vtk_file.isEqualString(float_str, "float")) {
854 for (Integer z = 0; z < nb_node_z; ++z) {
855 for (Integer y = 0; y < nb_node_y; ++y) {
856 for (Integer x = 0; x < nb_node_x; ++x) {
857 Real nx = vtk_file.getFloat();
858 Real ny = vtk_file.getFloat();
859 Real nz = vtk_file.getFloat();
860 Integer node_unique_id = x + y * nb_node_x + z * nb_node_xy;
861 coords[node_unique_id] = Real3(nx, ny, nz);
862 }
863 }
864 }
865 }
866 else if (vtk_file.isEqualString(float_str, "double")) {
867 for (Integer z = 0; z < nb_node_z; ++z) {
868 for (Integer y = 0; y < nb_node_y; ++y) {
869 for (Integer x = 0; x < nb_node_x; ++x) {
870 Real nx = vtk_file.getDouble();
871 Real ny = vtk_file.getDouble();
872 Real nz = vtk_file.getDouble();
873 Integer node_unique_id = x + y * nb_node_x + z * nb_node_xy;
874 coords[node_unique_id] = Real3(nx, ny, nz);
875 }
876 }
877 }
878 }
879 else {
880 throw IOException("_readStructuredGrid", "Invalid type name");
881 }
882
884 ENUMERATE_NODE (inode, mesh->allNodes()) {
885 Node node = *inode;
886 nodes_coord_var[inode] = coords[node.uniqueId().asInt32()];
887 }
888 }
889 }
890
891 // Créé les groupes de faces des côtés du parallélépipède
892 {
899
900 ENUMERATE_FACE (iface, mesh->allFaces()) {
901 Face face = *iface;
902 Integer face_local_id = face.localId();
903 bool is_xmin = true;
904 bool is_xmax = true;
905 bool is_ymin = true;
906 bool is_ymax = true;
907 bool is_zmin = true;
908 bool is_zmax = true;
909 for (Node node : face.nodes() ) {
910 Int64 node_unique_id = node.uniqueId().asInt64();
914 if (node_x != 0)
915 is_xmin = false;
916 if (node_x != (nb_node_x - 1))
917 is_xmax = false;
918 if (node_y != 0)
919 is_ymin = false;
920 if (node_y != (nb_node_y - 1))
921 is_ymax = false;
922 if (node_z != 0)
923 is_zmin = false;
924 if (node_z != (nb_node_z - 1))
925 is_zmax = false;
926 }
927 if (is_xmin)
928 xmin_surface_lid.add(face_local_id);
929 if (is_xmax)
930 xmax_surface_lid.add(face_local_id);
931 if (is_ymin)
932 ymin_surface_lid.add(face_local_id);
933 if (is_ymax)
934 ymax_surface_lid.add(face_local_id);
935 if (is_zmin)
936 zmin_surface_lid.add(face_local_id);
937 if (is_zmax)
938 zmax_surface_lid.add(face_local_id);
939 }
940 _createFaceGroup(mesh, "XMIN", xmin_surface_lid);
941 _createFaceGroup(mesh, "XMAX", xmax_surface_lid);
942 _createFaceGroup(mesh, "YMIN", ymin_surface_lid);
943 _createFaceGroup(mesh, "YMAX", ymax_surface_lid);
944 _createFaceGroup(mesh, "ZMIN", zmin_surface_lid);
945 _createFaceGroup(mesh, "ZMAX", zmax_surface_lid);
946 }
947
948 _readMetadata(mesh, vtk_file);
949
950 // Maintenant, regarde s'il existe des données associées aux fichier
952 if (r)
953 return r;
954
955 return false;
956}
957
958/*---------------------------------------------------------------------------*/
959/*---------------------------------------------------------------------------*/
960
970{
971 ARCANE_UNUSED(mesh);
972
973 const char* func_name = "VtkMeshIOService::_readNodesUnstructuredGrid()";
974 const char* buf = vtk_file.getNextLine();
975 std::istringstream iline(buf);
976 std::string points_str;
977 std::string data_type_str;
978 Integer nb_node = 0;
979
980 iline >> ws >> points_str >> ws >> nb_node >> ws >> data_type_str;
981
982 if (!iline)
983 throw IOException(func_name, "Syntax error while reading number of nodes");
984
985 vtk_file.checkString(points_str, "POINTS");
986
987 if (nb_node < 0)
988 throw IOException(A_FUNCINFO, String::format("Invalid number of nodes: n={0}", nb_node));
989
990 info() << " Info: " << nb_node;
991
992 // Lecture les coordonnées
993 node_coords.resize(nb_node);
994 {
995 if (vtk_file.isEqualString(data_type_str, "int")) {
996 for (Integer i = 0; i < nb_node; ++i) {
997 Real nx = vtk_file.getInt();
998 Real ny = vtk_file.getInt();
999 Real nz = vtk_file.getInt();
1000 node_coords[i] = Real3(nx, ny, nz);
1001 }
1002 }
1003 else if (vtk_file.isEqualString(data_type_str, "float")) {
1004 for (Integer i = 0; i < nb_node; ++i) {
1005 Real nx = vtk_file.getFloat();
1006 Real ny = vtk_file.getFloat();
1007 Real nz = vtk_file.getFloat();
1008 node_coords[i] = Real3(nx, ny, nz);
1009 }
1010 }
1011 else if (vtk_file.isEqualString(data_type_str, "double")) {
1012 for (Integer i = 0; i < nb_node; ++i) {
1013 Real nx = vtk_file.getDouble();
1014 Real ny = vtk_file.getDouble();
1015 Real nz = vtk_file.getDouble();
1016 node_coords[i] = Real3(nx, ny, nz);
1017 }
1018 }
1019 else {
1020 throw IOException(func_name, "Invalid type name");
1021 }
1022 }
1023 _readMetadata(mesh, vtk_file);
1024}
1025
1026/*---------------------------------------------------------------------------*/
1027/*---------------------------------------------------------------------------*/
1028
1042 Array<Integer>& cells_nb_node,
1043 Array<ItemTypeId>& cells_type,
1044 Array<Int64>& cells_connectivity)
1045{
1046 ARCANE_UNUSED(mesh);
1047
1048 const char* func_name = "VtkMeshIOService::_readCellsUnstructuredGrid()";
1049 const char* buf = vtk_file.getNextLine();
1050
1051 //String buftest = vtk_file.getCurrentLine(); // DEBUG
1052 //pinfo() << "Ligne lu : " << buftest.localstr();
1053
1054 std::istringstream iline(buf);
1055 std::string cells_str;
1056 Integer nb_cell = 0;
1057 Integer nb_cell_node = 0;
1058
1059 iline >> ws >> cells_str >> ws >> nb_cell >> ws >> nb_cell_node;
1060
1061 if (!iline)
1062 throw IOException(func_name, "Syntax error while reading cells");
1063
1064 vtk_file.checkString(cells_str, "CELLS");
1065
1066 if (nb_cell < 0 || nb_cell_node < 0) {
1067 throw IOException(A_FUNCINFO,
1068 String::format("Invalid dimensions: nb_cell={0} nb_cell_node={1}",
1069 nb_cell, nb_cell_node));
1070 }
1071
1072 cells_nb_node.resize(nb_cell);
1073 cells_type.resize(nb_cell);
1074 cells_connectivity.resize(nb_cell_node);
1075
1076 {
1077 Integer connectivity_index = 0;
1078 for (Integer i = 0; i < nb_cell; ++i) {
1079 Integer n = vtk_file.getInt();
1080 cells_nb_node[i] = n;
1081 for (Integer j = 0; j < n; ++j) {
1082 Integer id = vtk_file.getInt();
1083 cells_connectivity[connectivity_index] = id;
1085 }
1086 }
1087 }
1088
1089 _readMetadata(mesh, vtk_file);
1090
1091 // Lecture du type des mailles
1092 {
1093 buf = vtk_file.getNextLine();
1094 std::istringstream iline(buf);
1095 std::string cell_types_str;
1096 Integer nb_cell_type;
1097 iline >> ws >> cell_types_str >> ws >> nb_cell_type;
1098
1099 if (!iline) {
1100 throw IOException(func_name, "Syntax error while reading cell types");
1101 }
1102
1103 vtk_file.checkString(cell_types_str, "CELL_TYPES");
1104 if (nb_cell_type != nb_cell) {
1105 ARCANE_THROW(IOException, "Inconsistency in number of CELL_TYPES: v={0} nb_cell={1}",
1107 }
1108 }
1109
1110 for (Integer i = 0; i < nb_cell; ++i) {
1111 Integer vtk_ct = vtk_file.getInt();
1112 Int16 it = vtkToArcaneCellType(vtk_ct, cells_nb_node[i]);
1113 cells_type[i] = ItemTypeId{ it };
1114 }
1115 _readMetadata(mesh, vtk_file);
1116}
1117
1118/*---------------------------------------------------------------------------*/
1119/*---------------------------------------------------------------------------*/
1120
1130{
1131 ARCANE_UNUSED(mesh);
1132
1133 //const char* func_name = "VtkMeshIOService::_readMetadata()";
1134
1135 if (vtk_file.isEof())
1136 return false;
1137 String meta = vtk_file.getNextLine();
1138
1139 // METADATA ?
1140 if (!vtk_file.isEqualString(meta, "METADATA")) {
1141 // S'il n'y a pas de METADATA, on demande à ce que la ligne lue soit relue la prochaine fois.
1142 vtk_file.reReadSameLine();
1143 return false;
1144 }
1145
1146 // Tant qu'il n'y a pas de ligne vide, on lit.
1147 while (!vtk_file.isEmptyNextLine() && !vtk_file.isEof()) {
1148 }
1149 return false;
1150
1151 // // Si l'on a besoin de faire quelque chose avec les METADATA un jour, voilà un code non testé.
1152 // std::string trash;
1153 // Real trash_real;
1154 // const char* buf = vtk_file.getNextLine();
1155
1156 // // INFORMATION ou COMPONENT_NAMES
1157 // std::istringstream iline(buf);
1158
1159 // String name_str;
1160 // String data_type_str;
1161 // Integer nb_info = 0;
1162 // Integer size_vector = 0;
1163
1164 // iline >> ws >> name_str;
1165
1166 // if(vtk_file.isEqualString(name_str, "INFORMATION"))
1167 // {
1168 // iline >> ws >> nb_info;
1169 // for( Integer i = 0; i < nb_info; i++)
1170 // {
1171 // buf = vtk_file.getNextLine();
1172 // std::istringstream iline(buf);
1173
1174 // // NAME [key name] LOCATION [key location (e.g. class name)]
1175 // iline >> ws >> trash >> ws >> data_type_str >> ws >> trash >> ws >> name_str;
1176
1177 // buf = vtk_file.getNextLine();
1178 // std::istringstream iline(buf);
1179
1180 // if(vtk_file.isEqualString(data_type_str, "StringVector"))
1181 // {
1182 // iline >> ws >> trash >> ws >> size_vector;
1183 // for(Integer j = 0; j < size_vector; j++)
1184 // {
1185 // trash = vtk_file.getNextLine();
1186 // }
1187 // }
1188
1189 // // else if(vtk_file.isEqualString(data_type_str, "Double")
1190 // // || vtk_file.isEqualString(data_type_str, "IdType")
1191 // // || vtk_file.isEqualString(data_type_str, "Integer")
1192 // // || vtk_file.isEqualString(data_type_str, "UnsignedLong")
1193 // // || vtk_file.isEqualString(data_type_str, "String"))
1194 // // {
1195 // // iline >> ws >> trash >> ws >> trash_real;
1196 // // }
1197
1198 // // else if(vtk_file.isEqualString(data_type_str, "DoubleVector")
1199 // // || vtk_file.isEqualString(data_type_str, "IntegerVector")
1200 // // || vtk_file.isEqualString(data_type_str, "StringVector"))
1201 // // {
1202 // // iline >> ws >> trash >> ws >> size_vector;
1203 // // for(Integer j = 0; j < size_vector; j++)
1204 // // {
1205 // // iline >> ws >> trash_real;
1206 // // }
1207 // // }
1208 // }
1209 // }
1210
1211 // else if(vtk_file.isEqualString(name_str, "COMPONENT_NAMES"))
1212 // {
1213 // while(!vtk_file.isEmptyNextLine())
1214 // {
1215 // trash = vtk_file.getCurrentLine();
1216 // }
1217 // }
1218
1219 // else
1220 // {
1221 // throw IOException(func_name,"Syntax error after METADATA tag");
1222 // }
1223}
1224
1225/*---------------------------------------------------------------------------*/
1226/*---------------------------------------------------------------------------*/
1227
1238{
1239 Integer nb_node = 0;
1240 Integer nb_cell = 0;
1241 Integer nb_cell_node = 0;
1242 Int32 sid = mesh->parallelMng()->commRank();
1243 UniqueArray<Real3> node_coords;
1245
1246 // Si on utilise le partitionneur interne, seul le sous-domaine lit le maillage
1247 bool need_read = true;
1248
1250 need_read = (sid == 0);
1251
1252 std::array<Int64, 4> nb_cell_by_dimension = {};
1253 Int32 mesh_dimension = -1;
1254 ItemTypeMng* itm = mesh->itemTypeMng();
1256
1257 if (need_read) {
1258 // Lecture première partie du fichier (après header).
1259 _readNodesUnstructuredGrid(mesh, vtk_file, node_coords);
1260 debug() << "Lecture _readNodesUnstructuredGrid OK";
1261 nb_node = node_coords.size();
1262
1263 // Lecture des infos des mailles
1264 // Lecture de la connectivité
1265 UniqueArray<Integer> cells_nb_node;
1266 UniqueArray<Int64> cells_connectivity;
1267 UniqueArray<ItemTypeId> cells_type;
1268 _readCellsUnstructuredGrid(mesh, vtk_file, cells_nb_node, cells_type, cells_connectivity);
1269 debug() << "Lecture _readCellsUnstructuredGrid OK";
1270
1271 nb_cell = cells_nb_node.size();
1272 nb_cell_node = cells_connectivity.size();
1273 cells_local_id.resize(nb_cell);
1274
1275 // Création des mailles
1276 mesh_build_info.preAllocate(nb_cell, nb_cell_node);
1277
1278 {
1279 Int32 connectivity_index = 0;
1280 for (Integer i = 0; i < nb_cell; ++i) {
1281 Int32 current_cell_nb_node = cells_nb_node[i];
1282 Int64 cell_unique_id = i;
1283
1284 cells_local_id[i] = i;
1285
1286 Int16 cell_dim = itm->typeFromId(cells_type[i])->dimension();
1287 if (cell_dim >= 0 && cell_dim <= 3)
1289
1290 auto cell_nodes = cells_connectivity.subView(connectivity_index,current_cell_nb_node);
1291 mesh_build_info.addCell(cells_type[i],cell_unique_id, cell_nodes);
1293 }
1294 }
1295 // Vérifie qu'on n'a pas de mailles de différentes dimensions
1296 Int32 nb_different_dim = 0;
1297 for (Int32 i = 0; i < 4; ++i)
1298 if (nb_cell_by_dimension[i] != 0) {
1300 mesh_dimension = i;
1301 }
1302 if (nb_different_dim > 1)
1303 ARCANE_FATAL("The mesh contains cells of different dimension. nb0={0} nb1={1} nb2={2} nb3={3}",
1305 }
1306
1307 // Positionne la dimension du maillage.
1308 {
1310 wanted_dimension = mesh->parallelMng()->reduce(Parallel::ReduceMax, wanted_dimension);
1312 mesh_build_info.setMeshDimension(wanted_dimension);
1313 }
1314
1315 mesh_build_info.allocateMesh();
1316
1317 // Positionne les coordonnées
1318 {
1320 ENUMERATE_NODE (inode, mesh->allNodes()) {
1321 Node node = *inode;
1322 nodes_coord_var[inode] = node_coords[node.uniqueId().asInt32()];
1323 }
1324 }
1325
1326 // Maintenant, regarde s'il existe des données associées aux fichier
1328 debug() << "Lecture _readData OK";
1329
1330 return r;
1331}
1332
1333/*---------------------------------------------------------------------------*/
1334/*---------------------------------------------------------------------------*/
1335
1345_readFacesMesh(IMesh* mesh, const String& file_name, const String& dir_name,
1347{
1348 ARCANE_UNUSED(dir_name);
1349
1350 std::ifstream ifile(file_name.localstr(), std::ifstream::binary);
1351 if (!ifile) {
1352 info() << "No face descriptor file found '" << file_name << "'";
1353 return;
1354 }
1355
1357 const char* buf = 0;
1358
1359 // Lecture de la description
1360 String title = vtk_file.getNextLine();
1361 info() << "Titre du fichier VTK : " << title;
1362
1363 String format = vtk_file.getNextLine();
1364 if (VtkFile::isEqualString(format, "BINARY")) {
1365 vtk_file.setIsBinaryFile(true);
1366 }
1367
1368 eMeshType mesh_type = VTK_MT_Unknown;
1369 // Lecture du type de maillage
1370 // TODO: en parallèle, avec use_internal_partition vrai, seul le processeur 0
1371 // lit les données. Dans ce cas, inutile que les autres ouvre le fichier.
1372 {
1373 buf = vtk_file.getNextLine();
1374 std::istringstream mesh_type_line(buf);
1375 std::string dataset_str;
1376 std::string mesh_type_str;
1377 mesh_type_line >> ws >> dataset_str >> ws >> mesh_type_str;
1378 vtk_file.checkString(dataset_str, "DATASET");
1379
1380 if (VtkFile::isEqualString(mesh_type_str, "UNSTRUCTURED_GRID")) {
1381 mesh_type = VTK_MT_UnstructuredGrid;
1382 }
1383
1384 if (mesh_type == VTK_MT_Unknown) {
1385 error() << "Face descriptor file type must be 'UNSTRUCTURED_GRID' (format=" << mesh_type_str << "')";
1386 return;
1387 }
1388 }
1389
1390 {
1391 IParallelMng* pm = mesh->parallelMng();
1392 Integer nb_face = 0;
1393 Int32 sid = pm->commRank();
1394
1396
1397 // Si on utilise le partitionneur interne, seul le sous-domaine lit le maillage
1398 bool need_read = true;
1400 need_read = (sid == 0);
1401
1402 if (need_read) {
1403 {
1404 // Lit des noeuds, mais ne conserve pas leur coordonnées car cela n'est
1405 // pas nécessaire.
1406 UniqueArray<Real3> node_coords;
1407 _readNodesUnstructuredGrid(mesh, vtk_file, node_coords);
1408 //nb_node = node_coords.size();
1409 }
1410
1411 // Lecture des infos des faces
1412 // Lecture de la connectivité
1417 nb_face = faces_nb_node.size();
1418 //nb_face_node = faces_connectivity.size();
1419
1420 // Il faut à partir de la connectivité retrouver les localId() des faces
1421 faces_local_id.resize(nb_face);
1422 {
1423 IMeshUtilities* mu = mesh->utilities();
1424 mu->localIdsFromConnectivity(IK_Face, faces_nb_node, faces_connectivity, faces_local_id);
1425 }
1426 }
1427
1428 // Maintenant, regarde s'il existe des données associées aux fichiers
1430 }
1431}
1432
1433/*---------------------------------------------------------------------------*/
1434/*---------------------------------------------------------------------------*/
1435
1450{
1451 // Seul le sous-domain maitre lit les valeurs. Par contre, les autres
1452 // sous-domaines doivent connaitre la liste des variables et groupes créées.
1453 // Si une donnée porte le nom 'GROUP_*', on considère qu'il s'agit d'un
1454 // groupe
1455
1456 IParallelMng* pm = mesh->parallelMng();
1458
1459 Int32 sid = pm->commRank();
1460
1461 // Si pas de données, retourne immédiatement.
1462 {
1463 Byte has_data = 1;
1464 if ((sid == 0) && vtk_file.isEof())
1465 has_data = 0;
1466
1468 pm->broadcast(bb, 0);
1469 // Pas de data.
1470 if (!has_data)
1471 return false;
1472 }
1473
1475 created_infos_str() << "<?xml version='1.0' ?>\n";
1476 created_infos_str() << "<infos>";
1477
1478 Integer nb_cell_kind = mesh->nbItem(cell_kind);
1479 const char* buf = 0;
1480
1481 if (sid == 0) {
1482 bool reading_node = false;
1483 bool reading_cell = false;
1484 while (((buf = vtk_file.getNextLine()) != 0) && !vtk_file.isEof()) {
1485 debug() << "Read line";
1486 std::istringstream iline(buf);
1487 std::string data_str;
1488 iline >> data_str;
1489
1490 // Si l'on a un bloc "CELL_DATA".
1491 if (VtkFile::isEqualString(data_str, "CELL_DATA")) {
1492 Integer nb_item = 0;
1493 iline >> ws >> nb_item;
1494 reading_node = false;
1495 reading_cell = true;
1496 if (nb_item != nb_cell_kind)
1497 error() << "Size expected = " << nb_cell_kind << " found = " << nb_item;
1498 }
1499
1500 // Si l'on a un bloc "POINT_DATA".
1501 else if (VtkFile::isEqualString(data_str, "POINT_DATA")) {
1502 Integer nb_item = 0;
1503 iline >> ws >> nb_item;
1504 reading_node = true;
1505 reading_cell = false;
1506 if (nb_item != nb_node)
1507 error() << "Size expected = " << nb_node << " found = " << nb_item;
1508 }
1509
1510 // Si l'on a un bloc "FIELD".
1511 else if (VtkFile::isEqualString(data_str, "FIELD")) {
1512 std::string name_str;
1513 int nb_fields;
1514
1515 iline >> ws >> name_str >> ws >> nb_fields;
1516
1517 Integer nb_item = 0;
1518 std::string type_str;
1519 std::string s_name_str;
1520 int nb_component = 1;
1521 bool is_group = false;
1522
1523 for (Integer i = 0; i < nb_fields; i++) {
1524 buf = vtk_file.getNextLine();
1525 std::istringstream iline(buf);
1526 iline >> ws >> s_name_str >> ws >> nb_component >> ws >> nb_item >> ws >> type_str;
1527
1529 error() << "Size expected = " << nb_cell_kind << " found = " << nb_item;
1530
1532 error() << "Size expected = " << nb_node << " found = " << nb_item;
1533
1535 String cstr = name_str.substring(0, 6);
1536
1537 if (cstr == "GROUP_") {
1538 is_group = true;
1539 String new_name = name_str.substring(6);
1540 debug() << "** ** ** GROUP ! name=" << new_name;
1542 }
1543
1544 if (is_group) {
1545 if (!VtkFile::isEqualString(type_str, "int")) {
1546 error() << "Group type must be 'int', found=" << type_str;
1547 return true;
1548 }
1549
1550 if (reading_node) {
1551 created_infos_str() << "<node-group name='" << name_str << "'/>";
1553 }
1554
1555 if (reading_cell) {
1556 created_infos_str() << "<cell-group name='" << name_str << "'/>";
1558 }
1559 }
1560
1561 // TODO : Voir un exemple si possible.
1562 else {
1563 if (!VtkFile::isEqualString(type_str, "float") && !VtkFile::isEqualString(type_str, "double")) {
1564 error() << "Expecting 'float' or 'double' data type, found=" << type_str;
1565 return true;
1566 }
1567
1568 if (reading_node) {
1569 fatal() << "Unable to read POINT_DATA: feature not implemented";
1570 }
1571 if (reading_cell) {
1572 created_infos_str() << "<cell-variable name='" << name_str << "'/>";
1573
1574 if (cell_kind != IK_Cell)
1575 throw IOException("Unable to read face variables: feature not supported");
1576
1578 }
1579 }
1580 }
1581 }
1582
1583 else {
1584 // Lecture des valeurs (bloc "CELL_DATA" ou "POINT_DATA")
1585 if (reading_node || reading_cell) {
1586 std::string type_str;
1587 std::string s_name_str;
1588 //String name_str;
1589 bool is_group = false;
1590 int nb_component = 1;
1591
1592 iline >> ws >> s_name_str >> ws >> type_str >> ws >> nb_component;
1593 debug() << "** ** ** READNAME: name=" << s_name_str << " type=" << type_str;
1594
1596 String cstr = name_str.substring(0, 6);
1597
1598 if (cstr == "GROUP_") {
1599 is_group = true;
1600 String new_name = name_str.substring(6);
1601 info() << "** ** ** GROUP ! name=" << new_name;
1603 }
1604
1605 if (!VtkFile::isEqualString(data_str, "SCALARS")) {
1606 error() << "Expecting 'SCALARS' data type, found=" << data_str;
1607 return true;
1608 }
1609
1610 if (is_group) {
1611 if (!VtkFile::isEqualString(type_str, "int")) {
1612 error() << "Group type must be 'int', found=" << type_str;
1613 return true;
1614 }
1615
1616 // Pour lire LOOKUP_TABLE
1617 buf = vtk_file.getNextLine();
1618
1619 if (reading_node) {
1620 created_infos_str() << "<node-group name='" << name_str << "'/>";
1622 }
1623
1624 if (reading_cell) {
1625 created_infos_str() << "<cell-group name='" << name_str << "'/>";
1627 }
1628 }
1629 else {
1630 if (!VtkFile::isEqualString(type_str, "float") && !VtkFile::isEqualString(type_str, "double")) {
1631 error() << "Expecting 'float' or 'double' data type, found=" << type_str;
1632 return true;
1633 }
1634
1635 // Pour lire LOOKUP_TABLE
1636 /*buf = */ vtk_file.getNextLine();
1637 if (reading_node) {
1638 fatal() << "Unable to read POINT_DATA: feature not implemented";
1639 }
1640 if (reading_cell) {
1641 created_infos_str() << "<cell-variable name='" << name_str << "'/>";
1642
1643 if (cell_kind != IK_Cell)
1644 throw IOException("Unable to read face variables: feature not supported");
1645
1647 }
1648 }
1649 }
1650 else {
1651 error() << "Expecting value CELL_DATA or POINT_DATA, found='" << data_str << "'";
1652 return true;
1653 }
1654 }
1655 }
1656 }
1657 created_infos_str() << "</infos>";
1659 ByteUniqueArray bytes;
1660 if (sid == 0) {
1661 String str = created_infos_str.str();
1662 ByteConstArrayView bv = str.utf8();
1663 Integer len = bv.size();
1664 bytes.resize(len + 1);
1665 bytes.copy(bv);
1666 }
1667
1668 pm->broadcastMemoryBuffer(bytes, 0);
1669
1670 if (sid != 0) {
1671 String str = String::fromUtf8(bytes);
1672 info() << "FOUND STR=" << bytes.size() << " " << str;
1674 XmlNode doc_node = doc->documentNode();
1675
1676 // Lecture des variables
1677 {
1678 XmlNodeList vars = doc_node.documentElement().children("cell-variable");
1679 for (XmlNode xnode : vars) {
1680 String name = xnode.attrValue("name");
1681 info() << "Building variable: " << name;
1683 variable_mng->_internalApi()->addAutoDestroyVariable(var);
1684 }
1685 }
1686
1687 // Lecture des groupes de mailles
1688 {
1689 XmlNodeList vars = doc_node.documentElement().children("cell-group");
1691 for (XmlNode xnode : vars) {
1692 String name = xnode.attrValue("name");
1693 info() << "Building group: " << name;
1694 cell_family->createGroup(name);
1695 }
1696 }
1697
1698 // Lecture des groupes de noeuds
1699 {
1700 XmlNodeList vars = doc_node.documentElement().children("node-group");
1702 for (XmlNode xnode : vars) {
1703 String name = xnode.attrValue("name");
1704 info() << "Create node group: " << name;
1705 node_family->createGroup(name);
1706 }
1707 }
1708 }
1709 }
1710 return false;
1711}
1712
1713/*---------------------------------------------------------------------------*/
1714/*---------------------------------------------------------------------------*/
1715
1725{
1726 info() << "Building face group '" << name << "'"
1727 << " size=" << faces_lid.size();
1728
1729 mesh->faceFamily()->createGroup(name, faces_lid);
1730}
1731
1732/*---------------------------------------------------------------------------*/
1733/*---------------------------------------------------------------------------*/
1734
1745{
1746 //TODO Faire la conversion uniqueId() vers localId() correcte
1747 info() << "Reading values for variable: " << var_name << " n=" << nb_cell;
1748 auto* var = new VariableCellReal(VariableBuildInfo(mesh, var_name));
1749 mesh->variableMng()->_internalApi()->addAutoDestroyVariable(var);
1750 RealArrayView values(var->asArray());
1751 for (Integer i = 0; i < nb_cell; ++i) {
1752 Real v = vtk_file.getDouble();
1753 values[i] = v;
1754 }
1755 _readMetadata(mesh, vtk_file);
1756 info() << "Variable build finished: " << vtk_file.isEof();
1757}
1758
1759/*---------------------------------------------------------------------------*/
1760/*---------------------------------------------------------------------------*/
1761
1773_readItemGroup(IMesh* mesh, VtkFile& vtk_file, const String& name, Integer nb_item,
1775{
1777 info() << "Reading group info for group: " << name;
1778
1779 Int32UniqueArray ids;
1780 for (Integer i = 0; i < nb_item; ++i) {
1781 Integer v = vtk_file.getInt();
1782 if (v != 0)
1783 ids.add(local_id[i]);
1784 }
1785 info() << "Building group: " << name << " nb_element=" << ids.size();
1786
1787 item_family->createGroup(name, ids);
1788
1789 _readMetadata(mesh, vtk_file);
1790}
1791
1792/*---------------------------------------------------------------------------*/
1793/*---------------------------------------------------------------------------*/
1794
1804_readNodeGroup(IMesh* mesh, VtkFile& vtk_file, const String& name, Integer nb_item)
1805{
1807 info() << "Lecture infos groupes de noeuds pour le groupe: " << name;
1808
1809 Int32UniqueArray ids;
1810 for (Integer i = 0; i < nb_item; ++i) {
1811 Integer v = vtk_file.getInt();
1812 if (v != 0)
1813 ids.add(i);
1814 }
1815 info() << "Création groupe: " << name << " nb_element=" << ids.size();
1816
1817 item_family->createGroup(name, ids);
1818
1819 _readMetadata(mesh, vtk_file);
1820}
1821
1822/*---------------------------------------------------------------------------*/
1823/*---------------------------------------------------------------------------*/
1824
1825/*---------------------------------------------------------------------------*/
1826/*---------------------------------------------------------------------------*/
1827
1829: public BasicService
1830, public IMeshWriter
1831{
1832 public:
1833
1835 : BasicService(sbi)
1836 {}
1837
1838 public:
1839
1840 void build() override {}
1841
1842 public:
1843
1844 bool writeMeshToFile(IMesh* mesh, const String& file_name) override;
1845
1846 private:
1847
1849 void _saveGroups(IItemFamily* family, std::ostream& ofile);
1850};
1851
1852/*---------------------------------------------------------------------------*/
1853/*---------------------------------------------------------------------------*/
1854
1856
1857/*---------------------------------------------------------------------------*/
1858/*---------------------------------------------------------------------------*/
1875{
1878 _writeMeshToFile(mesh, fname + "faces.vtk", IK_Face);
1879 return false;
1880}
1881
1882/*---------------------------------------------------------------------------*/
1883/*---------------------------------------------------------------------------*/
1892{
1893 std::ofstream ofile(file_name.localstr());
1894 ofile.precision(FloatInfo<Real>::maxDigit());
1895 if (!ofile)
1896 throw IOException("VtkMeshIOService::writeMeshToFile(): Unable to open file");
1897 ofile << "# vtk DataFile Version 2.0\n";
1898 ofile << "Maillage Arcane\n";
1899 ofile << "ASCII\n";
1900 ofile << "DATASET UNSTRUCTURED_GRID\n";
1902 nodes_local_id_to_current.fill(NULL_ITEM_ID);
1903
1904 Integer nb_node = mesh->nbNode();
1906 Integer nb_cell_kind = cell_kind_family->nbItem();
1907
1908 // Sauve les nœuds
1909 {
1910 ofile << "POINTS " << nb_node << " double\n";
1911 VariableNodeReal3& coords(mesh->toPrimaryMesh()->nodesCoordinates());
1912 Integer node_index = 0;
1913 ENUMERATE_NODE (inode, mesh->allNodes()) {
1914 const Node& node = *inode;
1916 Real3 xyz = coords[inode];
1917 ofile << xyz.x << ' ' << xyz.y << ' ' << xyz.z << '\n';
1918 ++node_index;
1919 }
1920 }
1921
1922 // Sauve les mailles ou faces
1923 {
1925 ENUMERATE_ITEMWITHNODES(iitem, cell_kind_family->allItems())
1926 {
1927 nb_node_cell_kind += (*iitem).nbNode();
1928 }
1929 ofile << "CELLS " << nb_cell_kind << ' ' << nb_node_cell_kind << "\n";
1930 ENUMERATE_ITEMWITHNODES(iitem, cell_kind_family->allItems())
1931 {
1932 const ItemWithNodes& item = *iitem;
1933 Integer item_nb_node = item.nbNode();
1934 ofile << item_nb_node;
1935 for (NodeLocalId inode_lid : item.nodes() ) {
1937 }
1938 ofile << '\n';
1939 }
1940 // Le type doit être coherent avec celui de vtkCellType.h
1941 ofile << "CELL_TYPES " << nb_cell_kind << "\n";
1943 Int16 arcane_type = (*iitem).type();
1944 int type = arcaneToVtkCellType(arcane_type);
1945 ofile << type << '\n';
1946 }
1947 }
1948
1949 // Si on est dans le maillage des mailles, sauve les groupes de noeuds.
1950 if (cell_kind == IK_Cell) {
1951 ofile << "POINT_DATA " << nb_node << "\n";
1952 _saveGroups(mesh->itemFamily(IK_Node), ofile);
1953 }
1954
1955 // Sauve les groupes de mailles
1956 ofile << "CELL_DATA " << nb_cell_kind << "\n";
1957 _saveGroups(mesh->itemFamily(cell_kind), ofile);
1958}
1959
1960/*---------------------------------------------------------------------------*/
1961/*---------------------------------------------------------------------------*/
1962
1963void VtkLegacyMeshWriter::
1964_saveGroups(IItemFamily* family, std::ostream& ofile)
1965{
1966 info() << "Saving groups for family name=" << family->name();
1969 ItemGroup group = *igroup;
1970 // Inutile de sauver le groupe de toutes les entités
1971 if (group == family->allItems())
1972 continue;
1973 //HACK: a supprimer
1974 if (group.name() == "OuterFaces")
1975 continue;
1976 ofile << "SCALARS GROUP_" << group.name() << " int 1\n";
1977 ofile << "LOOKUP_TABLE default\n";
1978 in_group_list.fill('0');
1979 ENUMERATE_ITEM (iitem, group) {
1980 in_group_list[(*iitem).localId()] = '1';
1981 }
1982 ENUMERATE_ITEM (iitem, family->allItems()) {
1983 ofile << in_group_list[(*iitem).localId()] << '\n';
1984 }
1985 }
1986}
1987
1988/*---------------------------------------------------------------------------*/
1989/*---------------------------------------------------------------------------*/
1990
1991/*---------------------------------------------------------------------------*/
1992/*---------------------------------------------------------------------------*/
1993
1995: public AbstractService
1996, public IMeshReader
1997{
1998 public:
1999
2000 explicit VtkMeshReader(const ServiceBuildInfo& sbi)
2002 {}
2003
2004 public:
2005
2006 bool allowExtension(const String& str) override { return str == "vtk"; }
2008 const String& dir_name, bool use_internal_partition) override
2009
2010 {
2011 ARCANE_UNUSED(mesh_node);
2013 bool ret = vtk_service.readMesh(mesh, file_name, dir_name, use_internal_partition);
2014 if (ret)
2015 return RTError;
2016
2017 return RTOk;
2018 }
2019};
2020
2021/*---------------------------------------------------------------------------*/
2022/*---------------------------------------------------------------------------*/
2023
2025
2027 ServiceProperty("VtkLegacyMeshReader", ST_SubDomain),
2029
2030/*---------------------------------------------------------------------------*/
2031/*---------------------------------------------------------------------------*/
2032
2033/*---------------------------------------------------------------------------*/
2034/*---------------------------------------------------------------------------*/
2035
2037: public AbstractService
2038, public ICaseMeshReader
2039{
2040 public:
2041
2043 : public IMeshBuilder
2044 {
2045 public:
2046
2048 : m_trace_mng(tm)
2049 , m_read_info(read_info)
2050 {}
2051
2052 public:
2053
2055 {
2056 ARCANE_UNUSED(build_info);
2057 }
2059 {
2060 VtkMeshIOService vtk_service(m_trace_mng);
2061 String fname = m_read_info.fileName();
2062 m_trace_mng->info() << "VtkLegacy Reader (ICaseMeshReader) file_name=" << fname;
2063 bool ret = vtk_service.readMesh(pm, fname, m_read_info.directoryName(), m_read_info.isParallelRead());
2064 if (ret)
2065 ARCANE_FATAL("Can not read VTK File");
2066 }
2067
2068 private:
2069
2070 ITraceMng* m_trace_mng;
2071 CaseMeshReaderReadInfo m_read_info;
2072 };
2073
2074 public:
2075
2078 {}
2079
2080 public:
2081
2083 {
2084 IMeshBuilder* builder = nullptr;
2085 if (read_info.format() == "vtk")
2086 builder = new Builder(traceMng(), read_info);
2087 return makeRef(builder);
2088 }
2089};
2090
2091/*---------------------------------------------------------------------------*/
2092/*---------------------------------------------------------------------------*/
2093
2095 ServiceProperty("VtkLegacyCaseMeshReader", ST_SubDomain),
2097
2098/*---------------------------------------------------------------------------*/
2099/*---------------------------------------------------------------------------*/
2100
2101} // End namespace Arcane
2102
2103/*---------------------------------------------------------------------------*/
2104/*---------------------------------------------------------------------------*/
#define ARCANE_THROW(exception_class,...)
Macro pour envoyer une exception avec formattage.
#define ARCANE_FATAL(...)
Macro envoyant une exception FatalErrorException.
Types et macros pour itérer sur les entités du maillage.
#define ENUMERATE_FACE(name, group)
Enumérateur générique d'un groupe de faces.
#define ENUMERATE_(type, name, group)
Enumérateur générique d'un groupe d'entité
#define ENUMERATE_ITEM(name, group)
Enumérateur générique d'un groupe de noeuds.
#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.
#define ARCANE_SERVICE_INTERFACE(ainterface)
Macro pour déclarer une interface lors de l'enregistrement d'un service.
Classe de base d'un service.
Classe de base de service lié à un sous-domaine.
Informations nécessaires pour la lecture d'un fichier de maillage.
Face d'une maille.
Definition Item.h:932
Informations sur le type flottant.
Definition Limits.h:48
Interface du service de lecture du maillage à partir du jeu de données.
Interface d'une famille d'entités.
virtual ItemGroupCollection groups() const =0
Liste des groupes de cette famille.
virtual ItemGroup allItems() const =0
Groupe de toutes les entités.
virtual Int32 maxLocalId() const =0
virtual String name() const =0
Nom de la famille.
virtual IItemFamily * nodeFamily()=0
Retourne la famille des noeuds.
virtual FaceGroup allFaces()=0
Groupe de toutes les faces.
virtual IItemFamily * itemFamily(eItemKind ik)=0
Retourne la famille d'entité de type ik.
virtual Integer nbNode()=0
Nombre de noeuds du maillage.
virtual Integer nbItem(eItemKind ik)=0
Nombre d'éléments du genre ik.
virtual IItemFamily * faceFamily()=0
Retourne la famille des faces.
virtual NodeGroup allNodes()=0
Groupe de tous les noeuds.
Interface d'un service de création/lecture du maillage.
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
Interface d'une classe proposant des fonctions utilitaires sur maillage.
Interface d'un service d'écriture d'un maillage.
Definition IMeshWriter.h:36
virtual IParallelMng * parallelMng()=0
Gestionnaire de parallèlisme.
virtual IMeshUtilities * utilities()=0
Interface des fonctions utilitaires associée.
virtual ItemTypeMng * itemTypeMng() const =0
Gestionnaire de types d'entités associé
virtual IPrimaryMesh * toPrimaryMesh()=0
Retourne l'instance sous la forme d'un IPrimaryMesh.
virtual IVariableMng * variableMng() const =0
Gestionnaire de variable associé
Exception lorsqu'une erreur d'entrée/sortie est détectée.
Definition IOException.h:32
Interface du gestionnaire de parallélisme pour un sous-domaine.
virtual Int32 commRank() const =0
Rang de cette instance dans le communicateur.
virtual void broadcastMemoryBuffer(ByteArray &bytes, Int32 rank)=0
Effectue un broadcast d'une zone mémoire.
virtual VariableNodeReal3 & nodesCoordinates()=0
Coordonnées des noeuds.
virtual void allocateCells(Integer nb_cell, Int64ConstArrayView cells_infos, bool one_alloc=true)=0
Allocation d'un maillage.
virtual void endAllocate()=0
Indique une fin d'allocation de mailles.
virtual void setDimension(Integer dim)=0
Positionne la dimension du maillage (1D, 2D ou 3D).
Interface du gestionnaire de variables.
static IXmlDocumentHolder * loadFromBuffer(Span< const Byte > buffer, const String &name, ITraceMng *tm)
Charge un document XML.
Definition DomUtils.cc:426
Groupe d'entités de maillage.
Definition ItemGroup.h:49
const String & name() const
Nom du groupe.
Definition ItemGroup.h:76
Type d'une entité (Item).
Definition ItemTypeId.h:32
Gestionnaire des types d'entités de maillage.
Definition ItemTypeMng.h:66
Elément de maillage s'appuyant sur des noeuds (Edge,Face,Cell).
Definition Item.h:714
NodeConnectedListViewType nodes() const
Liste des noeuds de l'entité
Definition Item.h:771
Int32 nbNode() const
Nombre de noeuds de l'entité
Definition Item.h:765
constexpr Int32 localId() const
Identifiant local de l'entité dans le sous-domaine du processeur.
Definition Item.h:210
ItemUniqueId uniqueId() const
Identifiant unique sur tous les domaines.
Definition Item.h:216
Lecteur des fichiers de maillage via la bibliothèque LIMA.
Definition Lima.cc:149
Paramètres nécessaires à la construction d'un maillage.
Noeud d'un maillage.
Definition Item.h:564
Flot de sortie lié à une String.
Classe gérant un vecteur de réel de dimension 3.
Definition Real3.h:132
Structure contenant les informations pour créer un service.
Propriétés de création d'un service.
Informations pour allouer les entités d'un maillage non structuré.
Paramètres nécessaires à la construction d'une variable.
const char * getCurrentLine()
Permet de retourner la ligne présente dans le buffer.
bool isEmptyNextLine()
Permet de voir si la prochaine ligne est vide.
bool m_is_init
Y'a-t-il eu au moins une ligne lue.
std::istream * m_stream
Le stream.
static bool isEqualString(const String &current_value, const String &expected_value)
Permet de vérifier si expected_value == current_value.
double getDouble()
Permet de récupérer le double qui suit.
bool m_is_eof
Est-on à la fin du fichier.
char m_buf[BUFSIZE]
Le buffer contenant la ligne lue.
void getBinary(T &type)
Permet de récupérer le nombre binaire qui suit.
int getInt()
Permet de récupérer le int qui suit.
void checkString(const String &current_value, const String &expected_value)
Permet de vérifier si expected_value == current_value.
bool m_is_binary_file
Est-ce un fichier contenant des données en binaire.
bool m_need_reread_current_line
Doit-on relire la même ligne.
const char * getNextLine()
Permet de récupérer la prochaine ligne du fichier.
float getFloat()
Permet de récupérer le float qui suit.
void allocateMeshItems(IPrimaryMesh *pm) override
Alloue les entités du maillage géré par ce service.
void fillMeshBuildInfo(MeshBuildInfo &build_info) override
Remplit build_info avec les informations nécessaires pour créer le maillage.
Ref< IMeshBuilder > createBuilder(const CaseMeshReaderReadInfo &read_info) const override
Retourne un builder pour créer et lire le maillage dont les informations sont spécifiées dans read_in...
bool writeMeshToFile(IMesh *mesh, const String &file_name) override
Ecriture du maillage au Vtk.
void _writeMeshToFile(IMesh *mesh, const String &file_name, eItemKind cell_kind)
Ecrit le maillage au format Vtk.
void build() override
Construction de niveau build du service.
Lecteur des fichiers de maillage au format Vtk historique (legacy).
void _createFaceGroup(IMesh *mesh, const String &name, Int32ConstArrayView faces_lid)
Permet de créer un groupe de face de nom "name" et composé des faces ayant les ids inclus dans "faces...
bool _readMetadata(IMesh *mesh, VtkFile &vtk_file)
Lecture des metadata.
void _readCellVariable(IMesh *mesh, VtkFile &vtk_file, const String &name_str, Integer nb_cell)
Permet de créer une variable aux mailles à partir des infos du fichier vtk.
bool _readStructuredGrid(IPrimaryMesh *mesh, VtkFile &, bool use_internal_partition)
Permet de lire un fichier vtk contenant une STRUCTURED_GRID.
void _readNodesUnstructuredGrid(IMesh *mesh, VtkFile &vtk_file, Array< Real3 > &node_coords)
Lecture des noeuds et de leur coordonnées.
void _readCellsUnstructuredGrid(IMesh *mesh, VtkFile &vtk_file, Array< Integer > &cells_nb_node, Array< ItemTypeId > &cells_type, Array< Int64 > &cells_connectivity)
Lecture des mailles et de leur connectivité.
void _readItemGroup(IMesh *mesh, VtkFile &vtk_file, const String &name_str, Integer nb_item, eItemKind ik, ConstArrayView< Int32 > local_id)
Permet de créer un groupe d'item.
void _readFacesMesh(IMesh *mesh, const String &file_name, const String &dir_name, bool use_internal_partition)
Permet de lire le fichier truc.vtkfaces.vtk (s'il existe).
bool _readData(IMesh *mesh, VtkFile &vtk_file, bool use_internal_partition, eItemKind cell_kind, Int32ConstArrayView local_id, Integer nb_node)
Permet de lire les données complémentaires (POINT_DATA / CELL_DATA).
bool readMesh(IPrimaryMesh *mesh, const String &file_name, const String &dir_name, bool use_internal_partition)
Permet de débuter la lecture d'un fichier vtk.
void _readNodeGroup(IMesh *mesh, VtkFile &vtk_file, const String &name, Integer nb_item)
Permet de créer un groupe de node.
bool _readUnstructuredGrid(IPrimaryMesh *mesh, VtkFile &vtk_file, bool use_internal_partition)
Permet de lire un fichier vtk contenant une UNSTRUCTURED_GRID.
eReturnType readMeshFromFile(IPrimaryMesh *mesh, const XmlNode &mesh_node, const String &file_name, const String &dir_name, bool use_internal_partition) override
Lit un maillage à partir d'un fichier.
bool allowExtension(const String &str) override
Vérifie si le service supporte les fichiers avec l'extension str.
Liste de noeuds d'un arbre DOM.
Definition XmlNodeList.h:33
Noeud d'un arbre DOM.
Definition XmlNode.h:51
Integer size() const
Nombre d'éléments du vecteur.
Vue modifiable d'un tableau d'un type T.
void copy(Span< const T > rhs)
Copie les valeurs de rhs dans l'instance.
void resize(Int64 s)
Change le nombre d'éléments du tableau à s.
Vue constante d'un tableau de type T.
Interface du gestionnaire de traces.
virtual TraceMessage info()=0
Flot pour un message d'information.
Chaîne de caractères unicode.
ByteConstArrayView utf8() const
Retourne la conversion de l'instance dans l'encodage UTF-8.
Definition String.cc:275
const char * localstr() const
Retourne la conversion de l'instance dans l'encodage UTF-8.
Definition String.cc:227
ITraceMng * traceMng() const
Gestionnaire de trace.
TraceMessage error() const
Flot pour un message d'erreur.
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flot pour un message de debug.
TraceMessage info() const
Flot pour un message d'information.
TraceMessage fatal() const
Flot pour un message d'erreur fatale.
Vecteur 1D de données avec sémantique par valeur (style STL).
#define ARCANE_REGISTER_SERVICE(aclass, a_service_property,...)
Macro pour enregistrer un service.
MeshVariableScalarRefT< Cell, Real > VariableCellReal
Grandeur au centre des mailles de type réel.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
@ ST_SubDomain
Le service s'utilise au niveau du sous-domaine.
eItemKind
Genre d'entité de maillage.
@ IK_Node
Entité de maillage de genre noeud.
@ IK_Cell
Entité de maillage de genre maille.
@ IK_Face
Entité de maillage de genre face.