Arcane  v4.1.7.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-2026 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* VtkMeshIOService.cc (C) 2000-2026 */
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#include "arcane/utils/CheckedConvert.h"
27
28#include "arcane/core/BasicService.h"
29#include "arcane/core/FactoryService.h"
30#include "arcane/core/ICaseMeshReader.h"
31#include "arcane/core/IItemFamily.h"
32#include "arcane/core/IPrimaryMesh.h"
33#include "arcane/core/IMeshBuilder.h"
34#include "arcane/core/IMeshReader.h"
35#include "arcane/core/IMeshUtilities.h"
36#include "arcane/core/IMeshWriter.h"
37#include "arcane/core/IParallelMng.h"
38#include "arcane/core/IVariableAccessor.h"
39#include "arcane/core/IXmlDocumentHolder.h"
40#include "arcane/core/Item.h"
42#include "arcane/core/IVariableMng.h"
43#include "arcane/core/VariableTypes.h"
44#include "arcane/core/XmlNode.h"
45#include "arcane/core/XmlNodeList.h"
46#include "arcane/core/UnstructuredMeshAllocateBuildInfo.h"
47
48#include "arcane/core/internal/IVariableMngInternal.h"
49#include "arcane/core/internal/VtkCellTypes.h"
50
51/*---------------------------------------------------------------------------*/
52/*---------------------------------------------------------------------------*/
53
54namespace Arcane
55{
56
57/*---------------------------------------------------------------------------*/
58/*---------------------------------------------------------------------------*/
59
60class VtkFile;
61using namespace Arcane::VtkUtils;
62
63/*---------------------------------------------------------------------------*/
64/*---------------------------------------------------------------------------*/
102class VtkMeshIOService
103: public TraceAccessor
104{
105 public:
106
107 explicit VtkMeshIOService(ITraceMng* tm)
108 : TraceAccessor(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
141 bool readMesh(IPrimaryMesh* mesh, const String& file_name, const String& dir_name, bool use_internal_partition);
142
143 private:
144
145 bool _readStructuredGrid(IPrimaryMesh* mesh, VtkFile&, bool use_internal_partition);
146 bool _readUnstructuredGrid(IPrimaryMesh* mesh, VtkFile& vtk_file, bool use_internal_partition);
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,
149 eItemKind ik, ConstArrayView<Int32> local_id);
150 void _readNodeGroup(IMesh* mesh, VtkFile& vtk_file, const String& name, Integer nb_item);
151 void _createFaceGroup(IMesh* mesh, const String& name, Int32ConstArrayView faces_lid);
152 bool _readData(IMesh* mesh, VtkFile& vtk_file, bool use_internal_partition, eItemKind cell_kind,
153 Int32ConstArrayView local_id, Integer nb_node);
154 void _readNodesUnstructuredGrid(IMesh* mesh, VtkFile& vtk_file, Array<Real3>& node_coords);
156 Array<Int32>& cells_nb_node,
157 Array<ItemTypeId>& cells_type,
158 Array<Int64>& cells_connectivity);
159 void _readFacesMesh(IMesh* mesh, const String& file_name,
160 const String& dir_name, bool use_internal_partition);
161 bool _readMetadata(IMesh* mesh, VtkFile& vtk_file);
162
163 private:
164};
165
166/*---------------------------------------------------------------------------*/
167/*---------------------------------------------------------------------------*/
168
169class VtkFile
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
190 void checkString(const String& current_value, const String& expected_value);
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 ARCANE_THROW(IOException, "Unexpected EndOfFile");
272 }
273
274 if (!m_stream->good())
275 ARCANE_THROW(IOException, "Error when reading stream");
276
277 // Le getline s'arrete (par défaut) au char '\n' et ne l'inclus pas dans le buf
278 // mais le remplace par '\0'.
279 m_stream->getline(m_buf, sizeof(m_buf) - 1);
280
281 // Si on arrive au bout du fichier, on return true (pour dire oui, il y a une ligne vide,
282 // à l'appelant de gérer ça).
283 if (m_stream->eof()) {
284 m_is_eof = true;
285 return true;
286 }
287
288 // Sous Windows, une ligne vide commence par \r.
289 // getline remplace \n par \0, que ce soit sous Windows ou Linux.
290 if (m_buf[0] == '\r' || m_buf[0] == '\0') {
291 getNextLine();
292
293 // On demande à ce que le prochain appel à getNextLine renvoie la ligne
294 // qui vient tout juste d'être bufferisée.
296 return true;
297 }
298 else {
299 bool is_comment = true;
300
301 // On retire le commentaire, s'il y en a un, en remplaçant '#' par '\0'.
302 for (int i = 0; i < BUFSIZE && m_buf[i] != '\0'; ++i) {
303 if (!isspace(m_buf[i]) && m_buf[i] != '#' && is_comment) {
304 is_comment = false;
305 }
306 if (m_buf[i] == '#') {
307 m_buf[i] = '\0';
308 break;
309 }
310 }
311
312 // Si ce n'est pas un commentaire, on supprime juste le '\r' final (si windows).
313 if (!is_comment) {
314 // Supprime le '\r' final
315 for (int i = 0; i < BUFSIZE && m_buf[i] != '\0'; ++i) {
316 if (m_buf[i] == '\r') {
317 m_buf[i] = '\0';
318 break;
319 }
320 }
321 }
322
323 // Si c'était un commentaire, on recherche la prochaine ligne "valide"
324 // en appelant getNextLine.
325 else {
326 getNextLine();
327 }
328 }
330 return false;
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 ARCANE_THROW(IOException, "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 ARCANE_THROW(IOException, "Error when reading stream");
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 ARCANE_THROW(IOException, "Can not read '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 ARCANE_THROW(IOException, "Can not read '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 ARCANE_THROW(IOException, "Can not read '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.
490 Byte big_endian[sizeofT];
491 Byte little_endian[sizeofT];
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
519checkString(const String& current_value, const String& expected_value)
520{
521 String current_value_low = current_value.lower();
522 String expected_value_low = expected_value.lower();
523
524 if (current_value_low != expected_value_low) {
525 ARCANE_THROW(IOException, "Bad string. Expecting '{0}', buf found '{1}'",
526 expected_value, current_value);
527 }
528}
529
530/*---------------------------------------------------------------------------*/
531/*---------------------------------------------------------------------------*/
543checkString(const String& current_value, const String& expected_value1, const String& expected_value2)
544{
545 String current_value_low = current_value.lower();
546 String expected_value1_low = expected_value1.lower();
547 String expected_value2_low = expected_value2.lower();
548
549 if (current_value_low != expected_value1_low && current_value_low != expected_value2_low) {
550 ARCANE_THROW(IOException, "Bad string. Expecting '{0}' or '{1}, buf found '{2}'",
551 expected_value1, expected_value2, current_value);
552 }
553}
554
555/*---------------------------------------------------------------------------*/
556/*---------------------------------------------------------------------------*/
567isEqualString(const String& current_value, const String& expected_value)
568{
569 String current_value_low = current_value.lower();
570 String expected_value_low = expected_value.lower();
571 return (current_value_low == expected_value_low);
572}
573
574/*---------------------------------------------------------------------------*/
575/*---------------------------------------------------------------------------*/
576
577/*---------------------------------------------------------------------------*/
578/*---------------------------------------------------------------------------*/
579
590readMesh(IPrimaryMesh* mesh, const String& file_name, const String& dir_name, bool use_internal_partition)
591{
592 std::ifstream ifile(file_name.localstr(), std::ifstream::binary);
593
594 if (!ifile) {
595 error() << "Unable to read file '" << file_name << "'";
596 return true;
597 }
598
599 debug() << "Fichier ouvert : " << file_name.localstr();
600
601 VtkFile vtk_file(&ifile);
602 const char* buf = 0;
603
604 // Lecture de la description
605 // Lecture title.
606 String title = vtk_file.getNextLine();
607
608 info() << "Titre du fichier VTK : " << title.localstr();
609
610 // Lecture format.
611 String format = vtk_file.getNextLine();
612
613 debug() << "Format du fichier VTK : " << format.localstr();
614
615 if (VtkFile::isEqualString(format, "BINARY")) {
616 vtk_file.setIsBinaryFile(true);
617 }
618
619 eMeshType mesh_type = VTK_MT_Unknown;
620
621 // Lecture du type de maillage
622 // TODO: en parallèle, avec use_internal_partition vrai, seul le processeur 0
623 // lit les données. Dans ce cas, inutile que les autres ouvrent le fichier.
624 {
625 buf = vtk_file.getNextLine();
626
627 std::istringstream mesh_type_line(buf);
628 std::string dataset_str;
629 std::string mesh_type_str;
630
631 mesh_type_line >> ws >> dataset_str >> ws >> mesh_type_str;
632
633 vtk_file.checkString(dataset_str, "DATASET");
634
635 if (VtkFile::isEqualString(mesh_type_str, "STRUCTURED_GRID")) {
636 mesh_type = VTK_MT_StructuredGrid;
637 }
638
639 if (VtkFile::isEqualString(mesh_type_str, "UNSTRUCTURED_GRID")) {
640 mesh_type = VTK_MT_UnstructuredGrid;
641 }
642
643 if (mesh_type == VTK_MT_Unknown) {
644 error() << "Support exists only for 'STRUCTURED_GRID' and 'UNSTRUCTURED_GRID' formats (format=" << mesh_type_str << "')";
645 return true;
646 }
647 }
648 debug() << "Lecture en-tête OK";
649
650 bool ret = true;
651 switch (mesh_type) {
652 case VTK_MT_StructuredGrid:
653 ret = _readStructuredGrid(mesh, vtk_file, use_internal_partition);
654 break;
655
656 case VTK_MT_UnstructuredGrid:
657 ret = _readUnstructuredGrid(mesh, vtk_file, use_internal_partition);
658 debug() << "Lecture _readUnstructuredGrid OK";
659 if (!ret) {
660 // Tente de lire le fichier des faces s'il existe
661 _readFacesMesh(mesh, file_name + "faces.vtk", dir_name, use_internal_partition);
662 debug() << "Lecture _readFacesMesh OK";
663 }
664 break;
665
666 case VTK_MT_Unknown:
667 break;
668 }
669 /*while ( (buf=vtk_file.getNextLine()) != 0 ){
670 info() << " STR " << buf;
671 }*/
672
673 ifile.close();
674 return ret;
675}
676
677/*---------------------------------------------------------------------------*/
678/*---------------------------------------------------------------------------*/
679
689_readStructuredGrid(IPrimaryMesh* mesh, VtkFile& vtk_file, bool use_internal_partition)
690{
691 // Lecture du nombre de points: DIMENSIONS nx ny nz
692 const char* buf = nullptr;
693 Integer nb_node_x = 0;
694 Integer nb_node_y = 0;
695 Integer nb_node_z = 0;
696 {
697 buf = vtk_file.getNextLine();
698 std::istringstream iline(buf);
699 std::string dimension_str;
700 iline >> ws >> dimension_str >> ws >> nb_node_x >> ws >> nb_node_y >> ws >> nb_node_z;
701
702 if (!iline) {
703 error() << "Syntax error while reading grid dimensions";
704 return true;
705 }
706
707 vtk_file.checkString(dimension_str, "DIMENSIONS");
708 if (nb_node_x <= 1 || nb_node_y <= 1 || nb_node_z <= 1) {
709 error() << "Invalid dimensions: x=" << nb_node_x << " y=" << nb_node_y << " z=" << nb_node_z;
710 return true;
711 }
712 }
713 info() << " Infos: " << nb_node_x << " " << nb_node_y << " " << nb_node_z;
714 Integer nb_node = nb_node_x * nb_node_y * nb_node_z;
715
716 // Lecture du nombre de points: POINTS nb float
717 std::string float_str;
718 {
719 buf = vtk_file.getNextLine();
720 std::istringstream iline(buf);
721 std::string points_str;
722 Integer nb_node_read = 0;
723 iline >> ws >> points_str >> ws >> nb_node_read >> ws >> float_str;
724 if (!iline) {
725 error() << "Syntax error while reading grid dimensions";
726 return true;
727 }
728 vtk_file.checkString(points_str, "POINTS");
729 if (nb_node_read != nb_node) {
730 error() << "Number of invalid nodes: expected=" << nb_node << " found=" << nb_node_read;
731 return true;
732 }
733 }
734
735 Int32 rank = mesh->parallelMng()->commRank();
736
737 Integer nb_cell_x = nb_node_x - 1;
738 Integer nb_cell_y = nb_node_y - 1;
739 Integer nb_cell_z = nb_node_z - 1;
740
741 if (use_internal_partition && rank != 0) {
742 nb_node_x = 0;
743 nb_node_y = 0;
744 nb_node_z = 0;
745 nb_cell_x = 0;
746 nb_cell_y = 0;
747 nb_cell_z = 0;
748 }
749
750 const Integer nb_node_yz = nb_node_y * nb_node_z;
751 const Integer nb_node_xy = nb_node_x * nb_node_y;
752
753 Integer nb_cell = nb_cell_x * nb_cell_y * nb_cell_z;
754 UniqueArray<Int32> cells_local_id(nb_cell);
755
756 // Creation du maillage
757 {
758 UniqueArray<Integer> nodes_unique_id(nb_node);
759
760 info() << " NODE YZ = " << nb_node_yz;
761 // Création des noeuds
762 //Integer nb_node_local_id = 0;
763 {
764 Integer node_local_id = 0;
765 for (Integer x = 0; x < nb_node_x; ++x) {
766 for (Integer z = 0; z < nb_node_z; ++z) {
767 for (Integer y = 0; y < nb_node_y; ++y) {
768
769 Integer node_unique_id = y + (z)*nb_node_y + x * nb_node_y * nb_node_z;
770
771 nodes_unique_id[node_local_id] = node_unique_id;
772
773 ++node_local_id;
774 }
775 }
776 }
777 //nb_node_local_id = node_local_id;
778 //warning() << " NB NODE LOCAL ID=" << node_local_id;
779 }
780
781 // Création des mailles
782
783 // Infos pour la création des mailles
784 // par maille: 1 pour son unique id,
785 // 1 pour son type,
786 // 8 pour chaque noeud
787 UniqueArray<Int64> cells_infos(nb_cell * 10);
788
789 {
790 Integer cell_local_id = 0;
791 Integer cells_infos_index = 0;
792
793 // Normalement ne doit pas arriver car les valeurs de nb_node_x et
794 // nb_node_y sont testées lors de la lecture.
795 if (nb_node_xy == 0)
796 ARCANE_FATAL("Null value for nb_node_xy");
797
798 //Integer index = 0;
799 for (Integer z = 0; z < nb_cell_z; ++z) {
800 for (Integer y = 0; y < nb_cell_y; ++y) {
801 for (Integer x = 0; x < nb_cell_x; ++x) {
802 Integer current_cell_nb_node = 8;
803
804 //Integer cell_unique_id = y + (z)*nb_cell_y + x*nb_cell_y*nb_cell_z;
805 Int64 cell_unique_id = x + y * nb_cell_x + z * nb_cell_x * nb_cell_y;
806
807 cells_infos[cells_infos_index] = IT_Hexaedron8;
808 ++cells_infos_index;
809
810 cells_infos[cells_infos_index] = cell_unique_id;
811 ++cells_infos_index;
812
813 //Integer base_id = y + z*nb_node_y + x*nb_node_yz;
814 Integer base_id = x + y * nb_node_x + z * nb_node_xy;
815 cells_infos[cells_infos_index + 0] = nodes_unique_id[base_id];
816 cells_infos[cells_infos_index + 1] = nodes_unique_id[base_id + 1];
817 cells_infos[cells_infos_index + 2] = nodes_unique_id[base_id + nb_node_x + 1];
818 cells_infos[cells_infos_index + 3] = nodes_unique_id[base_id + nb_node_x + 0];
819 cells_infos[cells_infos_index + 4] = nodes_unique_id[base_id + nb_node_xy];
820 cells_infos[cells_infos_index + 5] = nodes_unique_id[base_id + nb_node_xy + 1];
821 cells_infos[cells_infos_index + 6] = nodes_unique_id[base_id + nb_node_xy + nb_node_x + 1];
822 cells_infos[cells_infos_index + 7] = nodes_unique_id[base_id + nb_node_xy + nb_node_x + 0];
823 cells_infos_index += current_cell_nb_node;
824 cells_local_id[cell_local_id] = cell_local_id;
825 ++cell_local_id;
826 }
827 }
828 }
829 }
830
831 mesh->setDimension(3);
832 mesh->allocateCells(nb_cell, cells_infos, false);
833 mesh->endAllocate();
834
835 // Positionne les coordonnées
836 {
837 UniqueArray<Real3> coords(nb_node);
838
839 if (vtk_file.isEqualString(float_str, "int")) {
840 for (Integer z = 0; z < nb_node_z; ++z) {
841 for (Integer y = 0; y < nb_node_y; ++y) {
842 for (Integer x = 0; x < nb_node_x; ++x) {
843 Real nx = vtk_file.getInt();
844 Real ny = vtk_file.getInt();
845 Real nz = vtk_file.getInt();
846 Integer node_unique_id = x + y * nb_node_x + z * nb_node_xy;
847 coords[node_unique_id] = Real3(nx, ny, nz);
848 }
849 }
850 }
851 }
852 else if (vtk_file.isEqualString(float_str, "float")) {
853 for (Integer z = 0; z < nb_node_z; ++z) {
854 for (Integer y = 0; y < nb_node_y; ++y) {
855 for (Integer x = 0; x < nb_node_x; ++x) {
856 Real nx = vtk_file.getFloat();
857 Real ny = vtk_file.getFloat();
858 Real nz = vtk_file.getFloat();
859 Integer node_unique_id = x + y * nb_node_x + z * nb_node_xy;
860 coords[node_unique_id] = Real3(nx, ny, nz);
861 }
862 }
863 }
864 }
865 else if (vtk_file.isEqualString(float_str, "double")) {
866 for (Integer z = 0; z < nb_node_z; ++z) {
867 for (Integer y = 0; y < nb_node_y; ++y) {
868 for (Integer x = 0; x < nb_node_x; ++x) {
869 Real nx = vtk_file.getDouble();
870 Real ny = vtk_file.getDouble();
871 Real nz = vtk_file.getDouble();
872 Integer node_unique_id = x + y * nb_node_x + z * nb_node_xy;
873 coords[node_unique_id] = Real3(nx, ny, nz);
874 }
875 }
876 }
877 }
878 else {
879 ARCANE_THROW(IOException, "Invalid type name");
880 }
881
882 VariableNodeReal3& nodes_coord_var(mesh->nodesCoordinates());
883 ENUMERATE_NODE (inode, mesh->allNodes()) {
884 Node node = *inode;
885 nodes_coord_var[inode] = coords[node.uniqueId().asInt32()];
886 }
887 }
888 }
889
890 // Créé les groupes de faces des côtés du parallélépipède
891 {
892 Int32UniqueArray xmin_surface_lid;
893 Int32UniqueArray xmax_surface_lid;
894 Int32UniqueArray ymin_surface_lid;
895 Int32UniqueArray ymax_surface_lid;
896 Int32UniqueArray zmin_surface_lid;
897 Int32UniqueArray zmax_surface_lid;
898
899 ENUMERATE_FACE (iface, mesh->allFaces()) {
900 Face face = *iface;
901 Integer face_local_id = face.localId();
902 bool is_xmin = true;
903 bool is_xmax = true;
904 bool is_ymin = true;
905 bool is_ymax = true;
906 bool is_zmin = true;
907 bool is_zmax = true;
908 for (Node node : face.nodes()) {
909 Int64 node_unique_id = node.uniqueId().asInt64();
910 Int64 node_z = node_unique_id / nb_node_xy;
911 Int64 node_y = (node_unique_id - node_z * nb_node_xy) / nb_node_x;
912 Int64 node_x = node_unique_id - node_z * nb_node_xy - node_y * nb_node_x;
913 if (node_x != 0)
914 is_xmin = false;
915 if (node_x != (nb_node_x - 1))
916 is_xmax = false;
917 if (node_y != 0)
918 is_ymin = false;
919 if (node_y != (nb_node_y - 1))
920 is_ymax = false;
921 if (node_z != 0)
922 is_zmin = false;
923 if (node_z != (nb_node_z - 1))
924 is_zmax = false;
925 }
926 if (is_xmin)
927 xmin_surface_lid.add(face_local_id);
928 if (is_xmax)
929 xmax_surface_lid.add(face_local_id);
930 if (is_ymin)
931 ymin_surface_lid.add(face_local_id);
932 if (is_ymax)
933 ymax_surface_lid.add(face_local_id);
934 if (is_zmin)
935 zmin_surface_lid.add(face_local_id);
936 if (is_zmax)
937 zmax_surface_lid.add(face_local_id);
938 }
939 _createFaceGroup(mesh, "XMIN", xmin_surface_lid);
940 _createFaceGroup(mesh, "XMAX", xmax_surface_lid);
941 _createFaceGroup(mesh, "YMIN", ymin_surface_lid);
942 _createFaceGroup(mesh, "YMAX", ymax_surface_lid);
943 _createFaceGroup(mesh, "ZMIN", zmin_surface_lid);
944 _createFaceGroup(mesh, "ZMAX", zmax_surface_lid);
945 }
946
947 _readMetadata(mesh, vtk_file);
948
949 // Maintenant, regarde s'il existe des données associées aux fichier
950 bool r = _readData(mesh, vtk_file, use_internal_partition, IK_Cell, cells_local_id, nb_node);
951 if (r)
952 return r;
953
954 return false;
955}
956
957/*---------------------------------------------------------------------------*/
958/*---------------------------------------------------------------------------*/
959
969{
970 ARCANE_UNUSED(mesh);
971
972 const char* buf = vtk_file.getNextLine();
973 std::istringstream iline(buf);
974 std::string points_str;
975 std::string data_type_str;
976 Integer nb_node = 0;
977
978 iline >> ws >> points_str >> ws >> nb_node >> ws >> data_type_str;
979
980 if (!iline)
981 ARCANE_THROW(IOException, "Syntax error while reading number of nodes");
982
983 vtk_file.checkString(points_str, "POINTS");
984
985 if (nb_node < 0)
986 ARCANE_THROW(IOException, "Invalid number of nodes: n={0}", nb_node);
987
988 info() << "VTK file : number of nodes = " << nb_node;
989
990 // Lecture les coordonnées
991 node_coords.resize(nb_node);
992 {
993 if (vtk_file.isEqualString(data_type_str, "int")) {
994 for (Integer i = 0; i < nb_node; ++i) {
995 Real nx = vtk_file.getInt();
996 Real ny = vtk_file.getInt();
997 Real nz = vtk_file.getInt();
998 node_coords[i] = Real3(nx, ny, nz);
999 }
1000 }
1001 else if (vtk_file.isEqualString(data_type_str, "float")) {
1002 for (Integer i = 0; i < nb_node; ++i) {
1003 Real nx = vtk_file.getFloat();
1004 Real ny = vtk_file.getFloat();
1005 Real nz = vtk_file.getFloat();
1006 node_coords[i] = Real3(nx, ny, nz);
1007 }
1008 }
1009 else if (vtk_file.isEqualString(data_type_str, "double")) {
1010 for (Integer i = 0; i < nb_node; ++i) {
1011 Real nx = vtk_file.getDouble();
1012 Real ny = vtk_file.getDouble();
1013 Real nz = vtk_file.getDouble();
1014 node_coords[i] = Real3(nx, ny, nz);
1015 }
1016 }
1017 else {
1018 ARCANE_THROW(IOException, "Invalid type name");
1019 }
1020 }
1021 _readMetadata(mesh, vtk_file);
1022}
1023
1024/*---------------------------------------------------------------------------*/
1025/*---------------------------------------------------------------------------*/
1026
1040 Array<Int32>& cells_nb_node,
1041 Array<ItemTypeId>& cells_type,
1042 Array<Int64>& cells_connectivity)
1043{
1044 ARCANE_UNUSED(mesh);
1045
1046 const char* buf = vtk_file.getNextLine();
1047
1048 std::istringstream iline(buf);
1049 std::string cells_str;
1050 Int32 i64_nb_cell = 0;
1051 Int32 i64_nb_cell_node = 0;
1052
1053 iline >> ws >> cells_str >> ws >> i64_nb_cell >> ws >> i64_nb_cell_node;
1054
1055 if (!iline)
1056 ARCANE_THROW(IOException, "Syntax error while reading cells");
1057
1058 vtk_file.checkString(cells_str, "CELLS");
1059
1060 info() << "VTK file : nb_cell = " << i64_nb_cell << " nb_cell_node=" << i64_nb_cell_node;
1061
1062 if (i64_nb_cell < 0 || i64_nb_cell_node < 0) {
1063 ARCANE_THROW(IOException, "Invalid dimensions: nb_cell={0} nb_cell_node={1}",
1064 i64_nb_cell, i64_nb_cell_node);
1065 }
1066
1067 Int32 nb_cell = CheckedConvert::toInt32(i64_nb_cell);
1068 Int32 nb_cell_node = CheckedConvert::toInt32(i64_nb_cell_node);
1069
1070 cells_nb_node.resize(nb_cell);
1071 cells_type.resize(nb_cell);
1072 cells_connectivity.resize(nb_cell_node);
1073
1074 {
1075 Integer connectivity_index = 0;
1076 for (Integer i = 0; i < nb_cell; ++i) {
1077 Integer n = vtk_file.getInt();
1078 cells_nb_node[i] = n;
1079 for (Integer j = 0; j < n; ++j) {
1080 Integer id = vtk_file.getInt();
1081 cells_connectivity[connectivity_index] = id;
1082 ++connectivity_index;
1083 }
1084 }
1085 }
1086
1087 _readMetadata(mesh, vtk_file);
1088
1089 // Lecture du type des mailles
1090 {
1091 buf = vtk_file.getNextLine();
1092 std::istringstream iline(buf);
1093 std::string cell_types_str;
1094 Integer nb_cell_type;
1095 iline >> ws >> cell_types_str >> ws >> nb_cell_type;
1096
1097 if (!iline) {
1098 ARCANE_THROW(IOException, "Syntax error while reading cell types");
1099 }
1100
1101 vtk_file.checkString(cell_types_str, "CELL_TYPES");
1102 if (nb_cell_type != nb_cell) {
1103 ARCANE_THROW(IOException, "Inconsistency in number of CELL_TYPES: v={0} nb_cell={1}",
1104 nb_cell_type, nb_cell);
1105 }
1106 }
1107
1108 for (Integer i = 0; i < nb_cell; ++i) {
1109 Integer vtk_ct = vtk_file.getInt();
1110 Int16 it = vtkToArcaneCellType(vtk_ct, cells_nb_node[i]);
1111 cells_type[i] = ItemTypeId{ it };
1112 }
1113 _readMetadata(mesh, vtk_file);
1114}
1115
1116/*---------------------------------------------------------------------------*/
1117/*---------------------------------------------------------------------------*/
1126_readMetadata(IMesh* mesh, VtkFile& vtk_file)
1127{
1128 ARCANE_UNUSED(mesh);
1129
1130 //const char* func_name = "VtkMeshIOService::_readMetadata()";
1131
1132 if (vtk_file.isEof())
1133 return false;
1134 String meta = vtk_file.getNextLine();
1135
1136 // METADATA ?
1137 if (!vtk_file.isEqualString(meta, "METADATA")) {
1138 // S'il n'y a pas de METADATA, on demande à ce que la ligne lue soit relue la prochaine fois.
1139 vtk_file.reReadSameLine();
1140 return false;
1141 }
1142
1143 // Tant qu'il n'y a pas de ligne vide, on lit.
1144 while (!vtk_file.isEmptyNextLine() && !vtk_file.isEof()) {
1145 }
1146 return false;
1147
1148 // // Si l'on a besoin de faire quelque chose avec les METADATA un jour, voilà un code non testé.
1149 // std::string trash;
1150 // Real trash_real;
1151 // const char* buf = vtk_file.getNextLine();
1152
1153 // // INFORMATION ou COMPONENT_NAMES
1154 // std::istringstream iline(buf);
1155
1156 // String name_str;
1157 // String data_type_str;
1158 // Integer nb_info = 0;
1159 // Integer size_vector = 0;
1160
1161 // iline >> ws >> name_str;
1162
1163 // if(vtk_file.isEqualString(name_str, "INFORMATION"))
1164 // {
1165 // iline >> ws >> nb_info;
1166 // for( Integer i = 0; i < nb_info; i++)
1167 // {
1168 // buf = vtk_file.getNextLine();
1169 // std::istringstream iline(buf);
1170
1171 // // NAME [key name] LOCATION [key location (e.g. class name)]
1172 // iline >> ws >> trash >> ws >> data_type_str >> ws >> trash >> ws >> name_str;
1173
1174 // buf = vtk_file.getNextLine();
1175 // std::istringstream iline(buf);
1176
1177 // if(vtk_file.isEqualString(data_type_str, "StringVector"))
1178 // {
1179 // iline >> ws >> trash >> ws >> size_vector;
1180 // for(Integer j = 0; j < size_vector; j++)
1181 // {
1182 // trash = vtk_file.getNextLine();
1183 // }
1184 // }
1185
1186 // // else if(vtk_file.isEqualString(data_type_str, "Double")
1187 // // || vtk_file.isEqualString(data_type_str, "IdType")
1188 // // || vtk_file.isEqualString(data_type_str, "Integer")
1189 // // || vtk_file.isEqualString(data_type_str, "UnsignedLong")
1190 // // || vtk_file.isEqualString(data_type_str, "String"))
1191 // // {
1192 // // iline >> ws >> trash >> ws >> trash_real;
1193 // // }
1194
1195 // // else if(vtk_file.isEqualString(data_type_str, "DoubleVector")
1196 // // || vtk_file.isEqualString(data_type_str, "IntegerVector")
1197 // // || vtk_file.isEqualString(data_type_str, "StringVector"))
1198 // // {
1199 // // iline >> ws >> trash >> ws >> size_vector;
1200 // // for(Integer j = 0; j < size_vector; j++)
1201 // // {
1202 // // iline >> ws >> trash_real;
1203 // // }
1204 // // }
1205 // }
1206 // }
1207
1208 // else if(vtk_file.isEqualString(name_str, "COMPONENT_NAMES"))
1209 // {
1210 // while(!vtk_file.isEmptyNextLine())
1211 // {
1212 // trash = vtk_file.getCurrentLine();
1213 // }
1214 // }
1215
1216 // else
1217 // {
1218 // throw IOException(func_name,"Syntax error after METADATA tag");
1219 // }
1220}
1221
1222/*---------------------------------------------------------------------------*/
1223/*---------------------------------------------------------------------------*/
1224
1234_readUnstructuredGrid(IPrimaryMesh* mesh, VtkFile& vtk_file, bool use_internal_partition)
1235{
1236 Integer nb_node = 0;
1237 Integer nb_cell = 0;
1238 Integer nb_cell_node = 0;
1239 Int32 sid = mesh->parallelMng()->commRank();
1240 UniqueArray<Real3> node_coords;
1241 UniqueArray<Int32> cells_local_id;
1242
1243 // Si on utilise le partitionneur interne, seul le sous-domaine lit le maillage
1244 bool need_read = true;
1245
1246 if (use_internal_partition)
1247 need_read = (sid == 0);
1248
1249 std::array<Int64, 4> nb_cell_by_dimension = {};
1250 Int32 mesh_dimension = -1;
1251 ItemTypeMng* itm = mesh->itemTypeMng();
1252 UnstructuredMeshAllocateBuildInfo mesh_build_info(mesh);
1253
1254 if (need_read) {
1255 // Lecture première partie du fichier (après header).
1256 _readNodesUnstructuredGrid(mesh, vtk_file, node_coords);
1257 debug() << "Lecture _readNodesUnstructuredGrid OK";
1258 nb_node = node_coords.size();
1259
1260 // Lecture des infos des mailles
1261 // Lecture de la connectivité
1262 UniqueArray<Int32> cells_nb_node;
1263 UniqueArray<Int64> cells_connectivity;
1264 UniqueArray<ItemTypeId> cells_type;
1265 _readCellsUnstructuredGrid(mesh, vtk_file, cells_nb_node, cells_type, cells_connectivity);
1266 debug() << "Lecture _readCellsUnstructuredGrid OK";
1267
1268 nb_cell = cells_nb_node.size();
1269 nb_cell_node = cells_connectivity.size();
1270 cells_local_id.resize(nb_cell);
1271
1272 // Création des mailles
1273 mesh_build_info.preAllocate(nb_cell, nb_cell_node);
1274
1275 {
1276 Int32 connectivity_index = 0;
1277 for (Integer i = 0; i < nb_cell; ++i) {
1278 Int32 current_cell_nb_node = cells_nb_node[i];
1279 Int64 cell_unique_id = i;
1280
1281 cells_local_id[i] = i;
1282
1283 Int16 cell_dim = itm->typeFromId(cells_type[i])->dimension();
1284 if (cell_dim >= 0 && cell_dim <= 3)
1285 ++nb_cell_by_dimension[cell_dim];
1286
1287 auto cell_nodes = cells_connectivity.subView(connectivity_index, current_cell_nb_node);
1288 mesh_build_info.addCell(cells_type[i], cell_unique_id, cell_nodes);
1289 connectivity_index += current_cell_nb_node;
1290 }
1291 }
1292 // Vérifie qu'on n'a pas de mailles de différentes dimensions
1293 Int32 nb_different_dim = 0;
1294 for (Int32 i = 0; i < 4; ++i)
1295 if (nb_cell_by_dimension[i] != 0) {
1296 ++nb_different_dim;
1297 mesh_dimension = i;
1298 }
1299 if (nb_different_dim > 1)
1300 ARCANE_FATAL("The mesh contains cells of different dimension. nb0={0} nb1={1} nb2={2} nb3={3}",
1301 nb_cell_by_dimension[0], nb_cell_by_dimension[1], nb_cell_by_dimension[2], nb_cell_by_dimension[3]);
1302 }
1303
1304 // Positionne la dimension du maillage.
1305 {
1306 Integer wanted_dimension = mesh_dimension;
1307 wanted_dimension = mesh->parallelMng()->reduce(Parallel::ReduceMax, wanted_dimension);
1308 //mesh->setDimension(wanted_dimension);
1309 mesh_build_info.setMeshDimension(wanted_dimension);
1310 }
1311
1312 mesh_build_info.allocateMesh();
1313
1314 // Positionne les coordonnées
1315 {
1316 VariableNodeReal3& nodes_coord_var(mesh->nodesCoordinates());
1317 ENUMERATE_NODE (inode, mesh->allNodes()) {
1318 Node node = *inode;
1319 nodes_coord_var[inode] = node_coords[node.uniqueId().asInt32()];
1320 }
1321 }
1322
1323 // Maintenant, regarde s'il existe des données associées aux fichier
1324 bool r = _readData(mesh, vtk_file, use_internal_partition, IK_Cell, cells_local_id, nb_node);
1325 debug() << "Lecture _readData OK";
1326
1327 return r;
1328}
1329
1330/*---------------------------------------------------------------------------*/
1331/*---------------------------------------------------------------------------*/
1332
1342_readFacesMesh(IMesh* mesh, const String& file_name, const String& dir_name,
1343 bool use_internal_partition)
1344{
1345 ARCANE_UNUSED(dir_name);
1346
1347 std::ifstream ifile(file_name.localstr(), std::ifstream::binary);
1348 if (!ifile) {
1349 info() << "No face descriptor file found '" << file_name << "'";
1350 return;
1351 }
1352
1353 VtkFile vtk_file(&ifile);
1354 const char* buf = 0;
1355
1356 // Lecture de la description
1357 String title = vtk_file.getNextLine();
1358 info() << "Reading VTK file '" << file_name << "'";
1359 info() << "Title of VTK file: " << title;
1360
1361 String format = vtk_file.getNextLine();
1362 if (VtkFile::isEqualString(format, "BINARY")) {
1363 vtk_file.setIsBinaryFile(true);
1364 }
1365
1366 eMeshType mesh_type = VTK_MT_Unknown;
1367 // Lecture du type de maillage
1368 // TODO: en parallèle, avec use_internal_partition vrai, seul le processeur 0
1369 // lit les données. Dans ce cas, inutile que les autres ouvre le fichier.
1370 {
1371 buf = vtk_file.getNextLine();
1372 std::istringstream mesh_type_line(buf);
1373 std::string dataset_str;
1374 std::string mesh_type_str;
1375 mesh_type_line >> ws >> dataset_str >> ws >> mesh_type_str;
1376 vtk_file.checkString(dataset_str, "DATASET");
1377
1378 if (VtkFile::isEqualString(mesh_type_str, "UNSTRUCTURED_GRID")) {
1379 mesh_type = VTK_MT_UnstructuredGrid;
1380 }
1381
1382 if (mesh_type == VTK_MT_Unknown) {
1383 error() << "Face descriptor file type must be 'UNSTRUCTURED_GRID' (format=" << mesh_type_str << "')";
1384 return;
1385 }
1386 }
1387
1388 {
1389 IParallelMng* pm = mesh->parallelMng();
1390 Integer nb_face = 0;
1391 Int32 sid = pm->commRank();
1392
1393 UniqueArray<Int32> faces_local_id;
1394
1395 // Si on utilise le partitionneur interne, seul le sous-domaine lit le maillage
1396 bool need_read = true;
1397 if (use_internal_partition)
1398 need_read = (sid == 0);
1399
1400 if (need_read) {
1401 {
1402 // Lit des noeuds, mais ne conserve pas leur coordonnées car cela n'est
1403 // pas nécessaire.
1404 UniqueArray<Real3> node_coords;
1405 _readNodesUnstructuredGrid(mesh, vtk_file, node_coords);
1406 //nb_node = node_coords.size();
1407 }
1408
1409 // Lecture des infos des faces
1410 // Lecture de la connectivité
1411 UniqueArray<Integer> faces_nb_node;
1412 UniqueArray<Int64> faces_connectivity;
1413 UniqueArray<ItemTypeId> faces_type;
1414 _readCellsUnstructuredGrid(mesh, vtk_file, faces_nb_node, faces_type, faces_connectivity);
1415 nb_face = faces_nb_node.size();
1416 //nb_face_node = faces_connectivity.size();
1417
1418 // Il faut à partir de la connectivité retrouver les localId() des faces
1419 faces_local_id.resize(nb_face);
1420 {
1421 IMeshUtilities* mu = mesh->utilities();
1422 mu->getFacesLocalIdFromConnectivity(faces_type, faces_connectivity, faces_local_id);
1423 }
1424 }
1425
1426 // Maintenant, regarde s'il existe des données associées aux fichiers
1427 _readData(mesh, vtk_file, use_internal_partition, IK_Face, faces_local_id, 0);
1428 }
1429}
1430
1431/*---------------------------------------------------------------------------*/
1432/*---------------------------------------------------------------------------*/
1433
1446_readData(IMesh* mesh, VtkFile& vtk_file, bool use_internal_partition,
1447 eItemKind cell_kind, Int32ConstArrayView local_id, Integer nb_node)
1448{
1449 // Seul le sous-domain maitre lit les valeurs. Par contre, les autres
1450 // sous-domaines doivent connaitre la liste des variables et groupes créées.
1451 // Si une donnée porte le nom 'GROUP_*', on considère qu'il s'agit d'un
1452 // groupe
1453
1454 IParallelMng* pm = mesh->parallelMng();
1455 IVariableMng* variable_mng = mesh->variableMng();
1456
1457 Int32 sid = pm->commRank();
1458
1459 // Si pas de données, retourne immédiatement.
1460 {
1461 Byte has_data = 1;
1462 if ((sid == 0) && vtk_file.isEof())
1463 has_data = 0;
1464
1465 ByteArrayView bb(1, &has_data);
1466 pm->broadcast(bb, 0);
1467 // Pas de data.
1468 if (!has_data)
1469 return false;
1470 }
1471
1472 OStringStream created_infos_str;
1473 created_infos_str() << "<?xml version='1.0' ?>\n";
1474 created_infos_str() << "<infos>";
1475
1476 Integer nb_cell_kind = mesh->nbItem(cell_kind);
1477 const char* buf = 0;
1478
1479 if (sid == 0) {
1480 bool reading_node = false;
1481 bool reading_cell = false;
1482 while (((buf = vtk_file.getNextLine()) != 0) && !vtk_file.isEof()) {
1483 debug() << "Read line";
1484 std::istringstream iline(buf);
1485 std::string data_str;
1486 iline >> data_str;
1487
1488 // Si l'on a un bloc "CELL_DATA".
1489 if (VtkFile::isEqualString(data_str, "CELL_DATA")) {
1490 Integer nb_item = 0;
1491 iline >> ws >> nb_item;
1492 reading_node = false;
1493 reading_cell = true;
1494 if (nb_item != nb_cell_kind)
1495 error() << "Size expected = " << nb_cell_kind << " found = " << nb_item;
1496 }
1497
1498 // Si l'on a un bloc "POINT_DATA".
1499 else if (VtkFile::isEqualString(data_str, "POINT_DATA")) {
1500 Integer nb_item = 0;
1501 iline >> ws >> nb_item;
1502 reading_node = true;
1503 reading_cell = false;
1504 if (nb_item != nb_node)
1505 error() << "Size expected = " << nb_node << " found = " << nb_item;
1506 }
1507
1508 // Si l'on a un bloc "FIELD".
1509 else if (VtkFile::isEqualString(data_str, "FIELD")) {
1510 std::string name_str;
1511 int nb_fields;
1512
1513 iline >> ws >> name_str >> ws >> nb_fields;
1514
1515 Integer nb_item = 0;
1516 std::string type_str;
1517 std::string s_name_str;
1518 int nb_component = 1;
1519 bool is_group = false;
1520
1521 for (Integer i = 0; i < nb_fields; i++) {
1522 buf = vtk_file.getNextLine();
1523 std::istringstream iline(buf);
1524 iline >> ws >> s_name_str >> ws >> nb_component >> ws >> nb_item >> ws >> type_str;
1525
1526 if (nb_item != nb_cell_kind && reading_cell && !reading_node)
1527 error() << "Size expected = " << nb_cell_kind << " found = " << nb_item;
1528
1529 if (nb_item != nb_node && !reading_cell && reading_node)
1530 error() << "Size expected = " << nb_node << " found = " << nb_item;
1531
1532 String name_str = s_name_str;
1533 String cstr = name_str.substring(0, 6);
1534
1535 if (cstr == "GROUP_") {
1536 is_group = true;
1537 String new_name = name_str.substring(6);
1538 debug() << "** ** ** GROUP ! name=" << new_name;
1539 name_str = new_name;
1540 }
1541
1542 if (is_group) {
1543 if (!VtkFile::isEqualString(type_str, "int")) {
1544 error() << "Group type must be 'int', found=" << type_str;
1545 return true;
1546 }
1547
1548 if (reading_node) {
1549 created_infos_str() << "<node-group name='" << name_str << "'/>";
1550 _readNodeGroup(mesh, vtk_file, name_str, nb_node);
1551 }
1552
1553 if (reading_cell) {
1554 created_infos_str() << "<cell-group name='" << name_str << "'/>";
1555 _readItemGroup(mesh, vtk_file, name_str, nb_cell_kind, cell_kind, local_id);
1556 }
1557 }
1558
1559 // TODO : Voir un exemple si possible.
1560 else {
1561 if (!VtkFile::isEqualString(type_str, "float") && !VtkFile::isEqualString(type_str, "double")) {
1562 error() << "Expecting 'float' or 'double' data type, found=" << type_str;
1563 return true;
1564 }
1565
1566 if (reading_node) {
1567 fatal() << "Unable to read POINT_DATA: feature not implemented";
1568 }
1569 if (reading_cell) {
1570 created_infos_str() << "<cell-variable name='" << name_str << "'/>";
1571
1572 if (cell_kind != IK_Cell)
1573 throw IOException("Unable to read face variables: feature not supported");
1574
1575 _readCellVariable(mesh, vtk_file, name_str, nb_cell_kind);
1576 }
1577 }
1578 }
1579 }
1580
1581 else {
1582 // Lecture des valeurs (bloc "CELL_DATA" ou "POINT_DATA")
1583 if (reading_node || reading_cell) {
1584 std::string type_str;
1585 std::string s_name_str;
1586 //String name_str;
1587 bool is_group = false;
1588 int nb_component = 1;
1589
1590 iline >> ws >> s_name_str >> ws >> type_str >> ws >> nb_component;
1591 debug() << "** ** ** READNAME: name=" << s_name_str << " type=" << type_str;
1592
1593 String name_str = s_name_str;
1594 String cstr = name_str.substring(0, 6);
1595
1596 if (cstr == "GROUP_") {
1597 is_group = true;
1598 String new_name = name_str.substring(6);
1599 info() << "** ** ** GROUP ! name=" << new_name;
1600 name_str = new_name;
1601 }
1602
1603 if (!VtkFile::isEqualString(data_str, "SCALARS")) {
1604 error() << "Expecting 'SCALARS' data type, found=" << data_str;
1605 return true;
1606 }
1607
1608 if (is_group) {
1609 if (!VtkFile::isEqualString(type_str, "int")) {
1610 error() << "Group type must be 'int', found=" << type_str;
1611 return true;
1612 }
1613
1614 // Pour lire LOOKUP_TABLE
1615 buf = vtk_file.getNextLine();
1616
1617 if (reading_node) {
1618 created_infos_str() << "<node-group name='" << name_str << "'/>";
1619 _readNodeGroup(mesh, vtk_file, name_str, nb_node);
1620 }
1621
1622 if (reading_cell) {
1623 created_infos_str() << "<cell-group name='" << name_str << "'/>";
1624 _readItemGroup(mesh, vtk_file, name_str, nb_cell_kind, cell_kind, local_id);
1625 }
1626 }
1627 else {
1628 if (!VtkFile::isEqualString(type_str, "float") && !VtkFile::isEqualString(type_str, "double")) {
1629 error() << "Expecting 'float' or 'double' data type, found=" << type_str;
1630 return true;
1631 }
1632
1633 // Pour lire LOOKUP_TABLE
1634 /*buf = */ vtk_file.getNextLine();
1635 if (reading_node) {
1636 fatal() << "Unable to read POINT_DATA: feature not implemented";
1637 }
1638 if (reading_cell) {
1639 created_infos_str() << "<cell-variable name='" << name_str << "'/>";
1640
1641 if (cell_kind != IK_Cell)
1642 throw IOException("Unable to read face variables: feature not supported");
1643
1644 _readCellVariable(mesh, vtk_file, name_str, nb_cell_kind);
1645 }
1646 }
1647 }
1648 else {
1649 error() << "Expecting value CELL_DATA or POINT_DATA, found='" << data_str << "'";
1650 return true;
1651 }
1652 }
1653 }
1654 }
1655 created_infos_str() << "</infos>";
1656 if (use_internal_partition) {
1657 ByteUniqueArray bytes;
1658 if (sid == 0) {
1659 String str = created_infos_str.str();
1660 ByteConstArrayView bv = str.utf8();
1661 Integer len = bv.size();
1662 bytes.resize(len + 1);
1663 bytes.copy(bv);
1664 }
1665
1666 pm->broadcastMemoryBuffer(bytes, 0);
1667
1668 if (sid != 0) {
1669 String str = String::fromUtf8(bytes);
1670 info() << "FOUND STR=" << bytes.size() << " " << str;
1672 XmlNode doc_node = doc->documentNode();
1673
1674 // Lecture des variables
1675 {
1676 XmlNodeList vars = doc_node.documentElement().children("cell-variable");
1677 for (XmlNode xnode : vars) {
1678 String name = xnode.attrValue("name");
1679 info() << "Building variable: " << name;
1681 variable_mng->_internalApi()->addAutoDestroyVariable(var);
1682 }
1683 }
1684
1685 // Lecture des groupes de mailles
1686 {
1687 XmlNodeList vars = doc_node.documentElement().children("cell-group");
1688 IItemFamily* cell_family = mesh->itemFamily(cell_kind);
1689 for (XmlNode xnode : vars) {
1690 String name = xnode.attrValue("name");
1691 info() << "Building group: " << name;
1692 cell_family->createGroup(name);
1693 }
1694 }
1695
1696 // Lecture des groupes de noeuds
1697 {
1698 XmlNodeList vars = doc_node.documentElement().children("node-group");
1699 IItemFamily* node_family = mesh->nodeFamily();
1700 for (XmlNode xnode : vars) {
1701 String name = xnode.attrValue("name");
1702 info() << "Create node group: " << name;
1703 node_family->createGroup(name);
1704 }
1705 }
1706 }
1707 }
1708 return false;
1709}
1710
1711/*---------------------------------------------------------------------------*/
1712/*---------------------------------------------------------------------------*/
1722_createFaceGroup(IMesh* mesh, const String& name, Int32ConstArrayView faces_lid)
1723{
1724 info() << "Building face group '" << name << "'"
1725 << " size=" << faces_lid.size();
1726
1727 mesh->faceFamily()->createGroup(name, faces_lid);
1728}
1729
1730/*---------------------------------------------------------------------------*/
1731/*---------------------------------------------------------------------------*/
1741_readCellVariable(IMesh* mesh, VtkFile& vtk_file, const String& var_name, Integer nb_cell)
1742{
1743 //TODO Faire la conversion uniqueId() vers localId() correcte
1744 info() << "Reading values for variable: " << var_name << " n=" << nb_cell;
1745 auto* var = new VariableCellReal(VariableBuildInfo(mesh, var_name));
1746 mesh->variableMng()->_internalApi()->addAutoDestroyVariable(var);
1747 RealArrayView values(var->asArray());
1748 for (Integer i = 0; i < nb_cell; ++i) {
1749 Real v = vtk_file.getDouble();
1750 values[i] = v;
1751 }
1752 _readMetadata(mesh, vtk_file);
1753 info() << "Variable build finished: " << vtk_file.isEof();
1754}
1755
1756/*---------------------------------------------------------------------------*/
1757/*---------------------------------------------------------------------------*/
1769_readItemGroup(IMesh* mesh, VtkFile& vtk_file, const String& name, Integer nb_item,
1770 eItemKind ik, Int32ConstArrayView local_id)
1771{
1772 IItemFamily* item_family = mesh->itemFamily(ik);
1773 info() << "Reading group info for group: " << name;
1774
1775 Int32UniqueArray ids;
1776 for (Integer i = 0; i < nb_item; ++i) {
1777 Integer v = vtk_file.getInt();
1778 if (v != 0)
1779 ids.add(local_id[i]);
1780 }
1781 info() << "Building group: " << name << " nb_element=" << ids.size();
1782
1783 item_family->createGroup(name, ids);
1784
1785 _readMetadata(mesh, vtk_file);
1786}
1787
1788/*---------------------------------------------------------------------------*/
1789/*---------------------------------------------------------------------------*/
1799_readNodeGroup(IMesh* mesh, VtkFile& vtk_file, const String& name, Integer nb_item)
1800{
1801 IItemFamily* item_family = mesh->itemFamily(IK_Node);
1802 info() << "Lecture infos groupes de noeuds pour le groupe: " << name;
1803
1804 Int32UniqueArray ids;
1805 for (Integer i = 0; i < nb_item; ++i) {
1806 Integer v = vtk_file.getInt();
1807 if (v != 0)
1808 ids.add(i);
1809 }
1810 info() << "Création groupe: " << name << " nb_element=" << ids.size();
1811
1812 item_family->createGroup(name, ids);
1813
1814 _readMetadata(mesh, vtk_file);
1815}
1816
1817/*---------------------------------------------------------------------------*/
1818/*---------------------------------------------------------------------------*/
1819
1820/*---------------------------------------------------------------------------*/
1821/*---------------------------------------------------------------------------*/
1822
1823class VtkLegacyMeshWriter
1824: public BasicService
1825, public IMeshWriter
1826{
1827 public:
1828
1829 explicit VtkLegacyMeshWriter(const ServiceBuildInfo& sbi)
1830 : BasicService(sbi)
1831 {}
1832
1833 public:
1834
1835 void build() override {}
1836
1837 public:
1838
1839 bool writeMeshToFile(IMesh* mesh, const String& file_name) override;
1840
1841 private:
1842
1843 void _writeMeshToFile(IMesh* mesh, const String& file_name, eItemKind cell_kind);
1844 void _saveGroups(IItemFamily* family, std::ostream& ofile);
1845};
1846
1847/*---------------------------------------------------------------------------*/
1848/*---------------------------------------------------------------------------*/
1849
1851 ServiceProperty("VtkLegacyMeshWriter", ST_SubDomain),
1853
1854/*---------------------------------------------------------------------------*/
1855/*---------------------------------------------------------------------------*/
1871writeMeshToFile(IMesh* mesh, const String& file_name)
1872{
1873 String fname = file_name;
1874 // Ajoute l'extension '.vtk' si elle n'y est pas.
1875 if (!fname.endsWith(".vtk"))
1876 fname = fname + ".vtk";
1877 _writeMeshToFile(mesh, fname, IK_Cell);
1878 _writeMeshToFile(mesh, fname + "faces.vtk", IK_Face);
1879 return false;
1880}
1881
1882/*---------------------------------------------------------------------------*/
1883/*---------------------------------------------------------------------------*/
1891_writeMeshToFile(IMesh* mesh, const String& file_name, eItemKind cell_kind)
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";
1901 UniqueArray<Integer> nodes_local_id_to_current(mesh->itemFamily(IK_Node)->maxLocalId());
1902 nodes_local_id_to_current.fill(NULL_ITEM_ID);
1903
1904 Integer nb_node = mesh->nbNode();
1905 IItemFamily* cell_kind_family = mesh->itemFamily(cell_kind);
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 Node node = *inode;
1915 nodes_local_id_to_current[node.localId()] = node_index;
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 {
1924 Integer nb_node_cell_kind = nb_cell_kind;
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 ItemWithNodes item = *iitem;
1933 Integer item_nb_node = item.nbNode();
1934 ofile << item_nb_node;
1935 for (NodeLocalId node_id : item.nodes()) {
1936 ofile << ' ' << nodes_local_id_to_current[node_id];
1937 }
1938 ofile << '\n';
1939 }
1940 // Le type doit être coherent avec celui de vtkCellType.h
1941 ofile << "CELL_TYPES " << nb_cell_kind << "\n";
1942 ENUMERATE_ (ItemWithNodes, iitem, cell_kind_family->allItems()) {
1943 int type = arcaneToVtkCellType(iitem->typeInfo());
1944 ofile << type << '\n';
1945 }
1946 }
1947
1948 // Si on est dans le maillage des mailles, sauve les groupes de noeuds.
1949 if (cell_kind == IK_Cell) {
1950 ofile << "POINT_DATA " << nb_node << "\n";
1951 _saveGroups(mesh->itemFamily(IK_Node), ofile);
1952 }
1953
1954 // Sauve les groupes de mailles
1955 ofile << "CELL_DATA " << nb_cell_kind << "\n";
1956 _saveGroups(mesh->itemFamily(cell_kind), ofile);
1957}
1958
1959/*---------------------------------------------------------------------------*/
1960/*---------------------------------------------------------------------------*/
1961
1962void VtkLegacyMeshWriter::
1963_saveGroups(IItemFamily* family, std::ostream& ofile)
1964{
1965 info() << "Saving groups for family name=" << family->name();
1966 UniqueArray<char> in_group_list(family->maxLocalId());
1967 for (ItemGroupCollection::Enumerator igroup(family->groups()); ++igroup;) {
1968 ItemGroup group = *igroup;
1969 // Inutile de sauver le groupe de toutes les entités
1970 if (group == family->allItems())
1971 continue;
1972 //HACK: a supprimer
1973 if (group.name() == "OuterFaces")
1974 continue;
1975 ofile << "SCALARS GROUP_" << group.name() << " int 1\n";
1976 ofile << "LOOKUP_TABLE default\n";
1977 in_group_list.fill('0');
1978 ENUMERATE_ITEM (iitem, group) {
1979 in_group_list[(*iitem).localId()] = '1';
1980 }
1981 ENUMERATE_ITEM (iitem, family->allItems()) {
1982 ofile << in_group_list[(*iitem).localId()] << '\n';
1983 }
1984 }
1985}
1986
1987/*---------------------------------------------------------------------------*/
1988/*---------------------------------------------------------------------------*/
1989
1990/*---------------------------------------------------------------------------*/
1991/*---------------------------------------------------------------------------*/
1992
1993class VtkMeshReader
1994: public AbstractService
1995, public IMeshReader
1996{
1997 public:
1998
1999 explicit VtkMeshReader(const ServiceBuildInfo& sbi)
2000 : AbstractService(sbi)
2001 {}
2002
2003 public:
2004
2005 bool allowExtension(const String& str) override { return str == "vtk"; }
2006 eReturnType readMeshFromFile(IPrimaryMesh* mesh, const XmlNode& mesh_node, const String& file_name,
2007 const String& dir_name, bool use_internal_partition) override
2008
2009 {
2010 ARCANE_UNUSED(mesh_node);
2011 VtkMeshIOService vtk_service(traceMng());
2012 bool ret = vtk_service.readMesh(mesh, file_name, dir_name, use_internal_partition);
2013 if (ret)
2014 return RTError;
2015
2016 return RTOk;
2017 }
2018};
2019
2020/*---------------------------------------------------------------------------*/
2021/*---------------------------------------------------------------------------*/
2022
2024
2026 ServiceProperty("VtkLegacyMeshReader", ST_SubDomain),
2028
2029/*---------------------------------------------------------------------------*/
2030/*---------------------------------------------------------------------------*/
2031
2032/*---------------------------------------------------------------------------*/
2033/*---------------------------------------------------------------------------*/
2034
2035class VtkLegacyCaseMeshReader
2036: public AbstractService
2037, public ICaseMeshReader
2038{
2039 public:
2040
2041 class Builder
2042 : public IMeshBuilder
2043 {
2044 public:
2045
2046 explicit Builder(ITraceMng* tm, const CaseMeshReaderReadInfo& read_info)
2047 : m_trace_mng(tm)
2048 , m_read_info(read_info)
2049 {}
2050
2051 public:
2052
2053 void fillMeshBuildInfo(MeshBuildInfo& build_info) override
2054 {
2055 ARCANE_UNUSED(build_info);
2056 }
2058 {
2059 VtkMeshIOService vtk_service(m_trace_mng);
2060 String fname = m_read_info.fileName();
2061 m_trace_mng->info() << "VtkLegacy Reader (ICaseMeshReader) file_name=" << fname;
2062 bool ret = vtk_service.readMesh(pm, fname, m_read_info.directoryName(), m_read_info.isParallelRead());
2063 if (ret)
2064 ARCANE_FATAL("Can not read VTK File");
2065 }
2066
2067 private:
2068
2069 ITraceMng* m_trace_mng;
2070 CaseMeshReaderReadInfo m_read_info;
2071 };
2072
2073 public:
2074
2075 explicit VtkLegacyCaseMeshReader(const ServiceBuildInfo& sbi)
2076 : AbstractService(sbi)
2077 {}
2078
2079 public:
2080
2082 {
2083 IMeshBuilder* builder = nullptr;
2084 if (read_info.format() == "vtk")
2085 builder = new Builder(traceMng(), read_info);
2086 return makeRef(builder);
2087 }
2088};
2089
2090/*---------------------------------------------------------------------------*/
2091/*---------------------------------------------------------------------------*/
2092
2094 ServiceProperty("VtkLegacyCaseMeshReader", ST_SubDomain),
2096
2097/*---------------------------------------------------------------------------*/
2098/*---------------------------------------------------------------------------*/
2099
2100} // End namespace Arcane
2101
2102/*---------------------------------------------------------------------------*/
2103/*---------------------------------------------------------------------------*/
#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.
Integer size() const
Nombre d'éléments du vecteur.
Classe de base d'un service.
AbstractService(const ServiceBuildInfo &)
Constructeur à partir d'un ServiceBuildInfo.
Tableau d'items de types quelconques.
void fill(const DataType &data)
Remplissage du tableau.
void resize(Int64 s)
Change le nombre d'éléments du tableau à s.
ArrayView< T > subView(Int64 abegin, Integer asize)
Sous-vue à partir de l'élément abegin et contenant asize éléments.
void copy(Span< const T > rhs)
Copie les valeurs de rhs dans l'instance.
void add(ConstReferenceType val)
Ajoute l'élément val à la fin du tableau.
Informations nécessaires pour la lecture d'un fichier de maillage.
Vue constante d'un tableau de type T.
constexpr Integer size() const noexcept
Nombre d'éléments du tableau.
Face d'une maille.
Definition Item.h:964
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.
Definition IItemFamily.h:84
virtual ItemGroupCollection groups() const =0
Liste des groupes de cette famille.
virtual ItemGroup allItems() const =0
Groupe de toutes les entités.
virtual ItemGroup createGroup(const String &name, Int32ConstArrayView local_ids, bool do_override=false)=0
Créé un groupe d'entités de nom name contenant les entités local_ids.
virtual Int32 maxLocalId() const =0
virtual String name() const =0
Nom de la famille.
virtual Integer nbItem() const =0
Nombre d'entités.
Interface d'un service de création/lecture du maillage.
Interface du service gérant la lecture d'un maillage.
Definition IMeshReader.h:32
eReturnType
Types des codes de retour d'une lecture ou écriture.
Definition IMeshReader.h:37
@ RTError
Erreur lors de l'opération.
Definition IMeshReader.h:39
@ RTOk
Opération effectuée avec succès.
Definition IMeshReader.h:38
Interface d'une classe proposant des fonctions utilitaires sur maillage.
virtual void getFacesLocalIdFromConnectivity(ConstArrayView< ItemTypeId > items_type, ConstArrayView< Int64 > items_connectivity, ArrayView< Int32 > local_ids, bool allow_null=false)=0
Recherche les identifiants locaux des faces à partir de leur connectivité.
Interface d'un service d'écriture d'un maillage.
Definition IMeshWriter.h:37
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.
Interface du gestionnaire de traces.
virtual void addAutoDestroyVariable(VariableRef *var)=0
Ajoute la variable à la liste des variables qui sont conservées jusqu'à la fin de l'exécution.
Interface du gestionnaire de variables.
virtual IVariableMngInternal * _internalApi()=0
API interne à Arcane.
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
Int16 dimension() const
Dimension de l'élément (<0 si inconnu)
Gestionnaire des types d'entités d'un maillage.
Definition ItemTypeMng.h:65
ItemTypeInfo * typeFromId(Integer id) const
Type correspondant au numéro id.
Elément de maillage s'appuyant sur des noeuds (Edge,Face,Cell).
Definition Item.h:736
NodeConnectedListViewType nodes() const
Liste des noeuds de l'entité
Definition Item.h:794
Int32 nbNode() const
Nombre de noeuds de l'entité
Definition Item.h:788
constexpr Int32 localId() const
Identifiant local de l'entité dans le sous-domaine du processeur.
Definition Item.h:219
ItemUniqueId uniqueId() const
Identifiant unique sur tous les domaines.
Definition Item.h:225
Paramètres nécessaires à la construction d'un maillage.
Noeud d'un maillage.
Definition Item.h:582
Flot de sortie lié à une String.
Classe gérant un vecteur de réel de dimension 3.
Definition Real3.h:132
Référence à une instance.
Encapsulation d'un pointeur qui se détruit automatiquement.
Definition ScopedPtr.h:44
Structure contenant les informations pour créer un service.
Propriétés de création d'un service.
Chaîne de caractères unicode.
String lower() const
Transforme tous les caractères de la chaîne en minuscules.
Definition String.cc:479
const char * localstr() const
Retourne la conversion de l'instance dans l'encodage UTF-8.
Definition String.cc:228
ByteConstArrayView utf8() const
Retourne la conversion de l'instance dans l'encodage UTF-8.
Definition String.cc:276
bool endsWith(const String &s) const
Indique si la chaîne se termine par les caractères de s.
Definition String.cc:1084
String substring(Int64 pos) const
Sous-chaîne commençant à la position pos.
Definition String.cc:1115
TraceAccessor(ITraceMng *m)
Construit un accesseur via le gestionnaire de trace m.
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flot pour un message de debug.
TraceMessage fatal() const
Flot pour un message d'erreur fatale.
TraceMessage info() const
Flot pour un message d'information.
TraceMessage error() const
Flot pour un message d'erreur.
ITraceMng * traceMng() const
Gestionnaire de trace.
Vecteur 1D de données avec sémantique par valeur (style STL).
Informations pour allouer les entités d'un maillage non structuré.
void allocateMesh()
Alloue le maillage avec les mailles ajoutées lors de l'appel à addCell().
void setMeshDimension(Int32 v)
Positionne la dimension du maillage.
void preAllocate(Int32 nb_cell, Int64 nb_connectivity_node)
Pre-alloue la mémoire.
void addCell(ItemTypeId type_id, Int64 cell_uid, SmallSpan< const Int64 > nodes_uid)
Ajoute une maille au maillage.
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 _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 _readCellsUnstructuredGrid(IMesh *mesh, VtkFile &vtk_file, Array< Int32 > &cells_nb_node, Array< ItemTypeId > &cells_type, Array< Int64 > &cells_connectivity)
Lecture des mailles et de leur connectivité.
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:35
Noeud d'un arbre DOM.
Definition XmlNode.h:51
XmlNode documentElement() const
Retourne le noeud élément du document.
Definition XmlNode.cc:556
XmlNodeList children(const String &name) const
Ensemble des noeuds fils de ce noeud ayant pour nom name.
Definition XmlNode.cc:93
#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.
MeshVariableScalarRefT< Node, Real3 > VariableNodeReal3
Grandeur au noeud de type coordonnées.
@ ReduceMax
Maximum des valeurs.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
ArrayView< Byte > ByteArrayView
Equivalent C d'un tableau à une dimension de caractères.
Definition UtilsTypes.h:447
std::int64_t Int64
Type entier signé sur 64 bits.
Int32 Integer
Type représentant un entier.
ConstArrayView< Int32 > Int32ConstArrayView
Equivalent C d'un tableau à une dimension d'entiers 32 bits.
Definition UtilsTypes.h:482
@ ST_SubDomain
Le service s'utilise au niveau du sous-domaine.
UniqueArray< Byte > ByteUniqueArray
Tableau dynamique à une dimension de caractères.
Definition UtilsTypes.h:335
UniqueArray< Int32 > Int32UniqueArray
Tableau dynamique à une dimension d'entiers 32 bits.
Definition UtilsTypes.h:341
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.
std::int16_t Int16
Type entier signé sur 16 bits.
double Real
Type représentant un réel.
ConstArrayView< Byte > ByteConstArrayView
Equivalent C d'un tableau à une dimension de caractères.
Definition UtilsTypes.h:476
unsigned char Byte
Type d'un octet.
Definition BaseTypes.h:43
auto makeRef(InstanceType *t) -> Ref< InstanceType >
Créé une référence sur un pointeur.
ArrayView< Real > RealArrayView
Equivalent C d'un tableau à une dimension de réels.
Definition UtilsTypes.h:459
std::int32_t Int32
Type entier signé sur 32 bits.
Real y
deuxième composante du triplet
Definition Real3.h:36
Real z
troisième composante du triplet
Definition Real3.h:37
Real x
première composante du triplet
Definition Real3.h:35