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