Arcane  4.1.12.0
Developer documentation
Loading...
Searching...
No Matches
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/* Reading/Writing a mesh in legacy VTK format. */
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
259{
260 m_is_init = true;
261
262 // We want getNextLine to read a new line.
263 // (we set it to false in case this method is
264 // called multiple times in a row).
266
267 // If we reached the end of the file during the previous call to this method or
268 // getNextLine, we throw an error.
269 if (m_is_eof) {
270 ARCANE_THROW(IOException, "Unexpected EndOfFile");
271 }
272
273 if (!m_stream->good())
274 ARCANE_THROW(IOException, "Error when reading stream");
275
276 // getline stops (by default) at the char '\n' and does not include it in the buf
277 // but replaces it with '\0'.
278 m_stream->getline(m_buf, sizeof(m_buf) - 1);
279
280 // If we reach the end of the file, we return true (to indicate yes, there is an empty line,
281 // for the caller to handle).
282 if (m_stream->eof()) {
283 m_is_eof = true;
284 return true;
285 }
286
287 // On Windows, an empty line starts with \r.
288 // getline replaces \n with \0, whether on Windows or Linux.
289 if (m_buf[0] == '\r' || m_buf[0] == '\0') {
290 getNextLine();
291
292 // We ask that the next call to getNextLine returns the line
293 // that was just buffered.
295 return true;
296 }
297 else {
298 bool is_comment = true;
299
300 // We remove the comment, if there is one, by replacing '#' with '\0'.
301 for (int i = 0; i < BUFSIZE && m_buf[i] != '\0'; ++i) {
302 if (!isspace(m_buf[i]) && m_buf[i] != '#' && is_comment) {
303 is_comment = false;
304 }
305 if (m_buf[i] == '#') {
306 m_buf[i] = '\0';
307 break;
308 }
309 }
310
311 // If it is not a comment, we just remove the final '\r' (if windows).
312 if (!is_comment) {
313 // Remove the final '\r'
314 for (int i = 0; i < BUFSIZE && m_buf[i] != '\0'; ++i) {
315 if (m_buf[i] == '\r') {
316 m_buf[i] = '\0';
317 break;
318 }
319 }
320 }
321
322 // If it was a comment, we search for the next "valid" line
323 // by calling getNextLine.
324 else {
325 getNextLine();
326 }
327 }
329 return false;
330}
331
332/*---------------------------------------------------------------------------*/
333/*---------------------------------------------------------------------------*/
334
340const char* VtkFile::
342{
343 m_is_init = true;
344
345 // We return the current buffer, if it has not been used.
348 return getCurrentLine();
349 }
350
351 // If we reached the end of the file during the previous call to this method or
352 // isEmptyNextLine, we throw an error.
353 if (m_is_eof) {
354 ARCANE_THROW(IOException, "Unexpected EndOfFile");
355 }
356
357 while (m_stream->good()) {
358 // getline stops (by default) at the char '\n' and does not include it in the buf but replaces it with '\0'.
359 m_stream->getline(m_buf, sizeof(m_buf) - 1);
360
361 // If we reach the end of the file, we return the buffer with \0 at the beginning (it is up to the caller to call
362 // isEof() to know if the file is finished or not).
363 if (m_stream->eof()) {
364 m_is_eof = true;
365 m_buf[0] = '\0';
366 return m_buf;
367 }
368
369 bool is_comment = true;
370
371 // On Windows, an empty line starts with \r.
372 // getline replaces \n with \0, whether on Windows or Linux.
373 if (m_buf[0] == '\0' || m_buf[0] == '\r')
374 continue;
375
376 // We remove the comment, if there is one, by replacing '#' with '\0'.
377 for (int i = 0; i < BUFSIZE && m_buf[i] != '\0'; ++i) {
378 if (!isspace(m_buf[i]) && m_buf[i] != '#' && is_comment) {
379 is_comment = false;
380 }
381 if (m_buf[i] == '#') {
382 m_buf[i] = '\0';
383 break;
384 }
385 }
386
387 // If it is not a comment, we just remove the final '\r' (if windows).
388 if (!is_comment) {
389 // Remove the final '\r'
390 for (int i = 0; i < BUFSIZE && m_buf[i] != '\0'; ++i) {
391 if (m_buf[i] == '\r') {
392 m_buf[i] = '\0';
393 break;
394 }
395 }
396 return m_buf;
397 }
398 }
399 ARCANE_THROW(IOException, "Error when reading stream");
400}
401
402/*---------------------------------------------------------------------------*/
403/*---------------------------------------------------------------------------*/
404
411getFloat()
412{
413 float v = 0.;
414 if (m_is_binary_file) {
415 getBinary(v);
416 return v;
417 }
418 (*m_stream) >> ws >> v;
419
420 if (m_stream->good())
421 return v;
422
423 ARCANE_THROW(IOException, "Can not read 'Float'");
424}
425
426/*---------------------------------------------------------------------------*/
427/*---------------------------------------------------------------------------*/
428
435getDouble()
436{
437 double v = 0.;
438 if (m_is_binary_file) {
439 getBinary(v);
440 return v;
441 }
442 (*m_stream) >> ws >> v;
443
444 if (m_stream->good())
445 return v;
446
447 ARCANE_THROW(IOException, "Can not read 'Double'");
448}
449
450/*---------------------------------------------------------------------------*/
451/*---------------------------------------------------------------------------*/
452
459getInt()
460{
461 int v = 0;
462 if (m_is_binary_file) {
463 getBinary(v);
464 return v;
465 }
466 (*m_stream) >> ws >> v;
467
468 if (m_stream->good())
469 return v;
470
471 ARCANE_THROW(IOException, "Can not read 'int'");
472}
473
474/*---------------------------------------------------------------------------*/
475/*---------------------------------------------------------------------------*/
476
482template <class T>
484getBinary(T& value)
485{
486 constexpr size_t sizeofT = sizeof(T);
487
488 // The VTK file is in big endian and current CPUs are in little endian.
489 Byte big_endian[sizeofT];
490 Byte little_endian[sizeofT];
491
492 // We read the next 'sizeofT' bytes and put them into big_endian.
493 m_stream->read((char*)big_endian, sizeofT);
494
495 // We transform big_endian into little_endian.
496 for (size_t i = 0; i < sizeofT; i++) {
497 little_endian[sizeofT - 1 - i] = big_endian[i];
498 }
499
500 // We 'cast' the byte array into type 'T'.
501 T* conv = new (little_endian) T;
502 value = *conv;
503}
504
505/*---------------------------------------------------------------------------*/
506/*---------------------------------------------------------------------------*/
507
518checkString(const String& current_value, const String& expected_value)
519{
520 String current_value_low = current_value.lower();
521 String expected_value_low = expected_value.lower();
522
523 if (current_value_low != expected_value_low) {
524 ARCANE_THROW(IOException, "Bad string. Expecting '{0}', buf found '{1}'",
525 expected_value, current_value);
526 }
527}
528
529/*---------------------------------------------------------------------------*/
530/*---------------------------------------------------------------------------*/
542checkString(const String& current_value, const String& expected_value1, const String& expected_value2)
543{
544 String current_value_low = current_value.lower();
545 String expected_value1_low = expected_value1.lower();
546 String expected_value2_low = expected_value2.lower();
547
548 if (current_value_low != expected_value1_low && current_value_low != expected_value2_low) {
549 ARCANE_THROW(IOException, "Bad string. Expecting '{0}' or '{1}, buf found '{2}'",
550 expected_value1, expected_value2, current_value);
551 }
552}
553
554/*---------------------------------------------------------------------------*/
555/*---------------------------------------------------------------------------*/
566isEqualString(const String& current_value, const String& expected_value)
567{
568 String current_value_low = current_value.lower();
569 String expected_value_low = expected_value.lower();
570 return (current_value_low == expected_value_low);
571}
572
573/*---------------------------------------------------------------------------*/
574/*---------------------------------------------------------------------------*/
575
586readMesh(IPrimaryMesh* mesh, const String& file_name, const String& dir_name, bool use_internal_partition)
587{
588 std::ifstream ifile(file_name.localstr(), std::ifstream::binary);
589
590 if (!ifile) {
591 error() << "Unable to read file '" << file_name << "'";
592 return true;
593 }
594
595 debug() << "Fichier ouvert : " << file_name.localstr();
596
597 VtkFile vtk_file(&ifile);
598 const char* buf = 0;
599
600 // Reading the description
601 // Reading title.
602 String title = vtk_file.getNextLine();
603
604 info() << "Titre du fichier VTK : " << title.localstr();
605
606 // Reading format.
607 String format = vtk_file.getNextLine();
608
609 debug() << "Format du fichier VTK : " << format.localstr();
610
611 if (VtkFile::isEqualString(format, "BINARY")) {
612 vtk_file.setIsBinaryFile(true);
613 }
614
615 eMeshType mesh_type = VTK_MT_Unknown;
616
617 // Reading the mesh type
618 // TODO: in parallel, with use_internal_partition true, only processor 0
619 // reads the data. In this case, it is unnecessary for others to open the file.
620 {
621 buf = vtk_file.getNextLine();
622
623 std::istringstream mesh_type_line(buf);
624 std::string dataset_str;
625 std::string mesh_type_str;
626
627 mesh_type_line >> ws >> dataset_str >> ws >> mesh_type_str;
628
629 vtk_file.checkString(dataset_str, "DATASET");
630
631 if (VtkFile::isEqualString(mesh_type_str, "STRUCTURED_GRID")) {
632 mesh_type = VTK_MT_StructuredGrid;
633 }
634
635 if (VtkFile::isEqualString(mesh_type_str, "UNSTRUCTURED_GRID")) {
636 mesh_type = VTK_MT_UnstructuredGrid;
637 }
638
639 if (mesh_type == VTK_MT_Unknown) {
640 error() << "Support exists only for 'STRUCTURED_GRID' and 'UNSTRUCTURED_GRID' formats (format=" << mesh_type_str << "')";
641 return true;
642 }
643 }
644 debug() << "Lecture en-tête OK";
645
646 bool ret = true;
647 switch (mesh_type) {
648 case VTK_MT_StructuredGrid:
649 ret = _readStructuredGrid(mesh, vtk_file, use_internal_partition);
650 break;
651
652 case VTK_MT_UnstructuredGrid:
653 ret = _readUnstructuredGrid(mesh, vtk_file, use_internal_partition);
654 debug() << "Lecture _readUnstructuredGrid OK";
655 if (!ret) {
656 // Tries to read the faces file if it exists
657 _readFacesMesh(mesh, file_name + "faces.vtk", dir_name, use_internal_partition);
658 debug() << "Lecture _readFacesMesh OK";
659 }
660 break;
661
662 case VTK_MT_Unknown:
663 break;
664 }
665 /*while ( (buf=vtk_file.getNextLine()) != 0 ){
666 info() << " STR " << buf;
667 }*/
668
669 ifile.close();
670 return ret;
671}
672
673/*---------------------------------------------------------------------------*/
674/*---------------------------------------------------------------------------*/
675
685_readStructuredGrid(IPrimaryMesh* mesh, VtkFile& vtk_file, bool use_internal_partition)
686{
687 // Reading the number of points: DIMENSIONS nx ny nz
688 const char* buf = nullptr;
689 Integer nb_node_x = 0;
690 Integer nb_node_y = 0;
691 Integer nb_node_z = 0;
692 {
693 buf = vtk_file.getNextLine();
694 std::istringstream iline(buf);
695 std::string dimension_str;
696 iline >> ws >> dimension_str >> ws >> nb_node_x >> ws >> nb_node_y >> ws >> nb_node_z;
697
698 if (!iline) {
699 error() << "Syntax error while reading grid dimensions";
700 return true;
701 }
702
703 vtk_file.checkString(dimension_str, "DIMENSIONS");
704 if (nb_node_x <= 1 || nb_node_y <= 1 || nb_node_z <= 1) {
705 error() << "Invalid dimensions: x=" << nb_node_x << " y=" << nb_node_y << " z=" << nb_node_z;
706 return true;
707 }
708 }
709 info() << " Infos: " << nb_node_x << " " << nb_node_y << " " << nb_node_z;
710 Integer nb_node = nb_node_x * nb_node_y * nb_node_z;
711
712 // Reading the number of points: POINTS nb float
713 std::string float_str;
714 {
715 buf = vtk_file.getNextLine();
716 std::istringstream iline(buf);
717 std::string points_str;
718 Integer nb_node_read = 0;
719 iline >> ws >> points_str >> ws >> nb_node_read >> ws >> float_str;
720 if (!iline) {
721 error() << "Syntax error while reading grid dimensions";
722 return true;
723 }
724 vtk_file.checkString(points_str, "POINTS");
725 if (nb_node_read != nb_node) {
726 error() << "Number of invalid nodes: expected=" << nb_node << " found=" << nb_node_read;
727 return true;
728 }
729 }
730
731 Int32 rank = mesh->parallelMng()->commRank();
732
733 Integer nb_cell_x = nb_node_x - 1;
734 Integer nb_cell_y = nb_node_y - 1;
735 Integer nb_cell_z = nb_node_z - 1;
736
737 if (use_internal_partition && rank != 0) {
738 nb_node_x = 0;
739 nb_node_y = 0;
740 nb_node_z = 0;
741 nb_cell_x = 0;
742 nb_cell_y = 0;
743 nb_cell_z = 0;
744 }
745
746 const Integer nb_node_yz = nb_node_y * nb_node_z;
747 const Integer nb_node_xy = nb_node_x * nb_node_y;
748
749 Integer nb_cell = nb_cell_x * nb_cell_y * nb_cell_z;
750 UniqueArray<Int32> cells_local_id(nb_cell);
751
752 // Mesh creation
753 {
754 UniqueArray<Integer> nodes_unique_id(nb_node);
755
756 info() << " NODE YZ = " << nb_node_yz;
757 // Node creation
758 //Integer nb_node_local_id = 0;
759 {
760 Integer node_local_id = 0;
761 for (Integer x = 0; x < nb_node_x; ++x) {
762 for (Integer z = 0; z < nb_node_z; ++z) {
763 for (Integer y = 0; y < nb_node_y; ++y) {
764
765 Integer node_unique_id = y + (z)*nb_node_y + x * nb_node_y * nb_node_z;
766
767 nodes_unique_id[node_local_id] = node_unique_id;
768
769 ++node_local_id;
770 }
771 }
772 }
773 //nb_node_local_id = node_local_id;
774 //warning() << " NB NODE LOCAL ID=" << node_local_id;
775 }
776
777 // Cell creation
778
779 // Info for cell creation
780 // per cell: 1 for its unique id,
781 // 1 for its type,
782 // 8 for each node
783 UniqueArray<Int64> cells_infos(nb_cell * 10);
784
785 {
786 Integer cell_local_id = 0;
787 Integer cells_infos_index = 0;
788
789 // Normally should not happen because the values of nb_node_x and
790 // nb_node_y are tested during reading.
791 if (nb_node_xy == 0)
792 ARCANE_FATAL("Null value for nb_node_xy");
793
794 //Integer index = 0;
795 for (Integer z = 0; z < nb_cell_z; ++z) {
796 for (Integer y = 0; y < nb_cell_y; ++y) {
797 for (Integer x = 0; x < nb_cell_x; ++x) {
798 Integer current_cell_nb_node = 8;
799
800 //Integer cell_unique_id = y + (z)*nb_cell_y + x*nb_cell_y*nb_cell_z;
801 Int64 cell_unique_id = x + y * nb_cell_x + z * nb_cell_x * nb_cell_y;
802
803 cells_infos[cells_infos_index] = IT_Hexaedron8;
804 ++cells_infos_index;
805
806 cells_infos[cells_infos_index] = cell_unique_id;
807 ++cells_infos_index;
808
809 //Integer base_id = y + z*nb_node_y + x*nb_node_yz;
810 Integer base_id = x + y * nb_node_x + z * nb_node_xy;
811 cells_infos[cells_infos_index + 0] = nodes_unique_id[base_id];
812 cells_infos[cells_infos_index + 1] = nodes_unique_id[base_id + 1];
813 cells_infos[cells_infos_index + 2] = nodes_unique_id[base_id + nb_node_x + 1];
814 cells_infos[cells_infos_index + 3] = nodes_unique_id[base_id + nb_node_x + 0];
815 cells_infos[cells_infos_index + 4] = nodes_unique_id[base_id + nb_node_xy];
816 cells_infos[cells_infos_index + 5] = nodes_unique_id[base_id + nb_node_xy + 1];
817 cells_infos[cells_infos_index + 6] = nodes_unique_id[base_id + nb_node_xy + nb_node_x + 1];
818 cells_infos[cells_infos_index + 7] = nodes_unique_id[base_id + nb_node_xy + nb_node_x + 0];
819 cells_infos_index += current_cell_nb_node;
820 cells_local_id[cell_local_id] = cell_local_id;
821 ++cell_local_id;
822 }
823 }
824 }
825 }
826
827 mesh->setDimension(3);
828 mesh->allocateCells(nb_cell, cells_infos, false);
829 mesh->endAllocate();
830
831 // Positioning the coordinates
832 {
833 UniqueArray<Real3> coords(nb_node);
834
835 if (vtk_file.isEqualString(float_str, "int")) {
836 for (Integer z = 0; z < nb_node_z; ++z) {
837 for (Integer y = 0; y < nb_node_y; ++y) {
838 for (Integer x = 0; x < nb_node_x; ++x) {
839 Real nx = vtk_file.getInt();
840 Real ny = vtk_file.getInt();
841 Real nz = vtk_file.getInt();
842 Integer node_unique_id = x + y * nb_node_x + z * nb_node_xy;
843 coords[node_unique_id] = Real3(nx, ny, nz);
844 }
845 }
846 }
847 }
848 else if (vtk_file.isEqualString(float_str, "float")) {
849 for (Integer z = 0; z < nb_node_z; ++z) {
850 for (Integer y = 0; y < nb_node_y; ++y) {
851 for (Integer x = 0; x < nb_node_x; ++x) {
852 Real nx = vtk_file.getFloat();
853 Real ny = vtk_file.getFloat();
854 Real nz = vtk_file.getFloat();
855 Integer node_unique_id = x + y * nb_node_x + z * nb_node_xy;
856 coords[node_unique_id] = Real3(nx, ny, nz);
857 }
858 }
859 }
860 }
861 else if (vtk_file.isEqualString(float_str, "double")) {
862 for (Integer z = 0; z < nb_node_z; ++z) {
863 for (Integer y = 0; y < nb_node_y; ++y) {
864 for (Integer x = 0; x < nb_node_x; ++x) {
865 Real nx = vtk_file.getDouble();
866 Real ny = vtk_file.getDouble();
867 Real nz = vtk_file.getDouble();
868 Integer node_unique_id = x + y * nb_node_x + z * nb_node_xy;
869 coords[node_unique_id] = Real3(nx, ny, nz);
870 }
871 }
872 }
873 }
874 else {
875 ARCANE_THROW(IOException, "Invalid type name");
876 }
877
878 VariableNodeReal3& nodes_coord_var(mesh->nodesCoordinates());
879 ENUMERATE_NODE (inode, mesh->allNodes()) {
880 Node node = *inode;
881 nodes_coord_var[inode] = coords[node.uniqueId().asInt32()];
882 }
883 }
884 }
885
886 // Created the face groups for the sides of the cuboid
887 {
888 Int32UniqueArray xmin_surface_lid;
889 Int32UniqueArray xmax_surface_lid;
890 Int32UniqueArray ymin_surface_lid;
891 Int32UniqueArray ymax_surface_lid;
892 Int32UniqueArray zmin_surface_lid;
893 Int32UniqueArray zmax_surface_lid;
894
895 ENUMERATE_FACE (iface, mesh->allFaces()) {
896 Face face = *iface;
897 Integer face_local_id = face.localId();
898 bool is_xmin = true;
899 bool is_xmax = true;
900 bool is_ymin = true;
901 bool is_ymax = true;
902 bool is_zmin = true;
903 bool is_zmax = true;
904 for (Node node : face.nodes()) {
905 Int64 node_unique_id = node.uniqueId().asInt64();
906 Int64 node_z = node_unique_id / nb_node_xy;
907 Int64 node_y = (node_unique_id - node_z * nb_node_xy) / nb_node_x;
908 Int64 node_x = node_unique_id - node_z * nb_node_xy - node_y * nb_node_x;
909 if (node_x != 0)
910 is_xmin = false;
911 if (node_x != (nb_node_x - 1))
912 is_xmax = false;
913 if (node_y != 0)
914 is_ymin = false;
915 if (node_y != (nb_node_y - 1))
916 is_ymax = false;
917 if (node_z != 0)
918 is_zmin = false;
919 if (node_z != (nb_node_z - 1))
920 is_zmax = false;
921 }
922 if (is_xmin)
923 xmin_surface_lid.add(face_local_id);
924 if (is_xmax)
925 xmax_surface_lid.add(face_local_id);
926 if (is_ymin)
927 ymin_surface_lid.add(face_local_id);
928 if (is_ymax)
929 ymax_surface_lid.add(face_local_id);
930 if (is_zmin)
931 zmin_surface_lid.add(face_local_id);
932 if (is_zmax)
933 zmax_surface_lid.add(face_local_id);
934 }
935 _createFaceGroup(mesh, "XMIN", xmin_surface_lid);
936 _createFaceGroup(mesh, "XMAX", xmax_surface_lid);
937 _createFaceGroup(mesh, "YMIN", ymin_surface_lid);
938 _createFaceGroup(mesh, "YMAX", ymax_surface_lid);
939 _createFaceGroup(mesh, "ZMIN", zmin_surface_lid);
940 _createFaceGroup(mesh, "ZMAX", zmax_surface_lid);
941 }
942
943 _readMetadata(mesh, vtk_file);
944
945 // Now, check if there is data associated with the file
946 bool r = _readData(mesh, vtk_file, use_internal_partition, IK_Cell, cells_local_id, nb_node);
947 if (r)
948 return r;
949
950 return false;
951}
952
953/*---------------------------------------------------------------------------*/
954/*---------------------------------------------------------------------------*/
955
965{
966 ARCANE_UNUSED(mesh);
967
968 const char* buf = vtk_file.getNextLine();
969 std::istringstream iline(buf);
970 std::string points_str;
971 std::string data_type_str;
972 Integer nb_node = 0;
973
974 iline >> ws >> points_str >> ws >> nb_node >> ws >> data_type_str;
975
976 if (!iline)
977 ARCANE_THROW(IOException, "Syntax error while reading number of nodes");
978
979 vtk_file.checkString(points_str, "POINTS");
980
981 if (nb_node < 0)
982 ARCANE_THROW(IOException, "Invalid number of nodes: n={0}", nb_node);
983
984 info() << "VTK file : number of nodes = " << nb_node;
985
986 // Read coordinates
987 node_coords.resize(nb_node);
988 {
989 if (vtk_file.isEqualString(data_type_str, "int")) {
990 for (Integer i = 0; i < nb_node; ++i) {
991 Real nx = vtk_file.getInt();
992 Real ny = vtk_file.getInt();
993 Real nz = vtk_file.getInt();
994 node_coords[i] = Real3(nx, ny, nz);
995 }
996 }
997 else if (vtk_file.isEqualString(data_type_str, "float")) {
998 for (Integer i = 0; i < nb_node; ++i) {
999 Real nx = vtk_file.getFloat();
1000 Real ny = vtk_file.getFloat();
1001 Real nz = vtk_file.getFloat();
1002 node_coords[i] = Real3(nx, ny, nz);
1003 }
1004 }
1005 else if (vtk_file.isEqualString(data_type_str, "double")) {
1006 for (Integer i = 0; i < nb_node; ++i) {
1007 Real nx = vtk_file.getDouble();
1008 Real ny = vtk_file.getDouble();
1009 Real nz = vtk_file.getDouble();
1010 node_coords[i] = Real3(nx, ny, nz);
1011 }
1012 }
1013 else {
1014 ARCANE_THROW(IOException, "Invalid type name");
1015 }
1016 }
1017 _readMetadata(mesh, vtk_file);
1018}
1019
1020/*---------------------------------------------------------------------------*/
1021/*---------------------------------------------------------------------------*/
1022
1036 Array<Int32>& cells_nb_node,
1037 Array<ItemTypeId>& cells_type,
1038 Array<Int64>& cells_connectivity)
1039{
1040 ARCANE_UNUSED(mesh);
1041
1042 const char* buf = vtk_file.getNextLine();
1043
1044 std::istringstream iline(buf);
1045 std::string cells_str;
1046 Int32 i64_nb_cell = 0;
1047 Int32 i64_nb_cell_node = 0;
1048
1049 iline >> ws >> cells_str >> ws >> i64_nb_cell >> ws >> i64_nb_cell_node;
1050
1051 if (!iline)
1052 ARCANE_THROW(IOException, "Syntax error while reading cells");
1053
1054 vtk_file.checkString(cells_str, "CELLS");
1055
1056 info() << "VTK file : nb_cell = " << i64_nb_cell << " nb_cell_node=" << i64_nb_cell_node;
1057
1058 if (i64_nb_cell < 0 || i64_nb_cell_node < 0) {
1059 ARCANE_THROW(IOException, "Invalid dimensions: nb_cell={0} nb_cell_node={1}",
1060 i64_nb_cell, i64_nb_cell_node);
1061 }
1062
1063 Int32 nb_cell = CheckedConvert::toInt32(i64_nb_cell);
1064 Int32 nb_cell_node = CheckedConvert::toInt32(i64_nb_cell_node);
1065
1066 cells_nb_node.resize(nb_cell);
1067 cells_type.resize(nb_cell);
1068 cells_connectivity.resize(nb_cell_node);
1069
1070 {
1071 Integer connectivity_index = 0;
1072 for (Integer i = 0; i < nb_cell; ++i) {
1073 Integer n = vtk_file.getInt();
1074 cells_nb_node[i] = n;
1075 for (Integer j = 0; j < n; ++j) {
1076 Integer id = vtk_file.getInt();
1077 cells_connectivity[connectivity_index] = id;
1078 ++connectivity_index;
1079 }
1080 }
1081 }
1082
1083 _readMetadata(mesh, vtk_file);
1084
1085 // Read cell types
1086 {
1087 buf = vtk_file.getNextLine();
1088 std::istringstream iline(buf);
1089 std::string cell_types_str;
1090 Integer nb_cell_type;
1091 iline >> ws >> cell_types_str >> ws >> nb_cell_type;
1092
1093 if (!iline) {
1094 ARCANE_THROW(IOException, "Syntax error while reading cell types");
1095 }
1096
1097 vtk_file.checkString(cell_types_str, "CELL_TYPES");
1098 if (nb_cell_type != nb_cell) {
1099 ARCANE_THROW(IOException, "Inconsistency in number of CELL_TYPES: v={0} nb_cell={1}",
1100 nb_cell_type, nb_cell);
1101 }
1102 }
1103
1104 for (Integer i = 0; i < nb_cell; ++i) {
1105 Integer vtk_ct = vtk_file.getInt();
1106 Int16 it = vtkToArcaneCellType(vtk_ct, cells_nb_node[i]);
1107 cells_type[i] = ItemTypeId{ it };
1108 }
1109 _readMetadata(mesh, vtk_file);
1110}
1111
1112/*---------------------------------------------------------------------------*/
1113/*---------------------------------------------------------------------------*/
1114
1123_readMetadata(IMesh* mesh, VtkFile& vtk_file)
1124{
1125 ARCANE_UNUSED(mesh);
1126
1127 //const char* func_name = "VtkMeshIOService::_readMetadata()";
1128
1129 if (vtk_file.isEof())
1130 return false;
1131 String meta = vtk_file.getNextLine();
1132
1133 // METADATA ?
1134 if (!vtk_file.isEqualString(meta, "METADATA")) {
1135 // If there is no METADATA, we ask that the line read be reread next time.
1136 vtk_file.reReadSameLine();
1137 return false;
1138 }
1139
1140 // As long as there is no empty line, we read.
1141 while (!vtk_file.isEmptyNextLine() && !vtk_file.isEof()) {
1142 }
1143 return false;
1144
1145 // // If we need to do something with METADATA someday, here is untested code.
1146 // std::string trash;
1147 // Real trash_real;
1148 // const char* buf = vtk_file.getNextLine();
1149
1150 // // INFORMATION or COMPONENT_NAMES
1151 // std::istringstream iline(buf);
1152
1153 // String name_str;
1154 // String data_type_str;
1155 // Integer nb_info = 0;
1156 // Integer size_vector = 0;
1157
1158 // iline >> ws >> name_str;
1159
1160 // if(vtk_file.isEqualString(name_str, "INFORMATION"))
1161 // {
1162 // iline >> ws >> nb_info;
1163 // for( Integer i = 0; i < nb_info; i++)
1164 // {
1165 // buf = vtk_file.getNextLine();
1166 // std::istringstream iline(buf);
1167
1168 // // NAME [key name] LOCATION [key location (e.g. class name)]
1169 // iline >> ws >> trash >> ws >> data_type_str >> ws >> trash >> ws >> name_str;
1170
1171 // buf = vtk_file.getNextLine();
1172 // std::istringstream iline(buf);
1173
1174 // if(vtk_file.isEqualString(data_type_str, "StringVector"))
1175 // {
1176 // iline >> ws >> trash >> ws >> size_vector;
1177 // for(Integer j = 0; j < size_vector; j++)
1178 // {
1179 // trash = vtk_file.getNextLine();
1180 // }
1181 // }
1182
1183 // // else if(vtk_file.isEqualString(data_type_str, "Double")
1184 // // || vtk_file.isEqualString(data_type_str, "IdType")
1185 // // || vtk_file.isEqualString(data_type_str, "Integer")
1186 // // || vtk_file.isEqualString(data_type_str, "UnsignedLong")
1187 // // || vtk_file.isEqualString(data_type_str, "String"))
1188 // // {
1189 // // iline >> ws >> trash >> ws >> trash_real;
1190 // // }
1191
1192 // // else if(vtk_file.isEqualString(data_type_str, "DoubleVector")
1193 // // || vtk_file.isEqualString(data_type_str, "IntegerVector")
1194 // // || vtk_file.isEqualString(data_type_str, "StringVector"))
1195 // // {
1196 // // iline >> ws >> trash >> ws >> size_vector;
1197 // // for(Integer j = 0; j < size_vector; j++)
1198 // // {
1199 // // iline >> ws >> trash_real;
1200 // // }
1201 // // }
1202 // }
1203 // }
1204
1205 // else if(vtk_file.isEqualString(name_str, "COMPONENT_NAMES"))
1206 // {
1207 // while(!vtk_file.isEmptyNextLine())
1208 // {
1209 // trash = vtk_file.getCurrentLine();
1210 // }
1211 // }
1212
1213 // else
1214 // {
1215 // throw IOException(func_name,"Syntax error after METADATA tag");
1216 // }
1217}
1218
1219/*---------------------------------------------------------------------------*/
1220/*---------------------------------------------------------------------------*/
1221
1231_readUnstructuredGrid(IPrimaryMesh* mesh, VtkFile& vtk_file, bool use_internal_partition)
1232{
1233 Integer nb_node = 0;
1234 Integer nb_cell = 0;
1235 Integer nb_cell_node = 0;
1236 Int32 sid = mesh->parallelMng()->commRank();
1237 UniqueArray<Real3> node_coords;
1238 UniqueArray<Int32> cells_local_id;
1239
1240 // If we use the internal partitioner, only the subdomain reads the mesh
1241 bool need_read = true;
1242
1243 if (use_internal_partition)
1244 need_read = (sid == 0);
1245
1246 std::array<Int64, 4> nb_cell_by_dimension = {};
1247 Int32 mesh_dimension = -1;
1248 ItemTypeMng* itm = mesh->itemTypeMng();
1249 UnstructuredMeshAllocateBuildInfo mesh_build_info(mesh);
1250
1251 if (need_read) {
1252 // Read the first part of the file (after header).
1253 _readNodesUnstructuredGrid(mesh, vtk_file, node_coords);
1254 debug() << "Lecture _readNodesUnstructuredGrid OK";
1255 nb_node = node_coords.size();
1256
1257 // Read cell info
1258 // Read connectivity
1259 UniqueArray<Int32> cells_nb_node;
1260 UniqueArray<Int64> cells_connectivity;
1261 UniqueArray<ItemTypeId> cells_type;
1262 _readCellsUnstructuredGrid(mesh, vtk_file, cells_nb_node, cells_type, cells_connectivity);
1263 debug() << "Reading _readCellsUnstructuredGrid OK";
1264
1265 nb_cell = cells_nb_node.size();
1266 nb_cell_node = cells_connectivity.size();
1267 cells_local_id.resize(nb_cell);
1268
1269 // Cell creation
1270 mesh_build_info.preAllocate(nb_cell, nb_cell_node);
1271
1272 {
1273 Int32 connectivity_index = 0;
1274 for (Integer i = 0; i < nb_cell; ++i) {
1275 Int32 current_cell_nb_node = cells_nb_node[i];
1276 Int64 cell_unique_id = i;
1277
1278 cells_local_id[i] = i;
1279
1280 Int16 cell_dim = itm->typeFromId(cells_type[i])->dimension();
1281 if (cell_dim >= 0 && cell_dim <= 3)
1282 ++nb_cell_by_dimension[cell_dim];
1283
1284 auto cell_nodes = cells_connectivity.subView(connectivity_index, current_cell_nb_node);
1285 mesh_build_info.addCell(cells_type[i], cell_unique_id, cell_nodes);
1286 connectivity_index += current_cell_nb_node;
1287 }
1288 }
1289 // Check that there are no cells of different dimensions
1290 Int32 nb_different_dim = 0;
1291 for (Int32 i = 0; i < 4; ++i)
1292 if (nb_cell_by_dimension[i] != 0) {
1293 ++nb_different_dim;
1294 mesh_dimension = i;
1295 }
1296 if (nb_different_dim > 1)
1297 ARCANE_FATAL("The mesh contains cells of different dimension. nb0={0} nb1={1} nb2={2} nb3={3}",
1298 nb_cell_by_dimension[0], nb_cell_by_dimension[1], nb_cell_by_dimension[2], nb_cell_by_dimension[3]);
1299 }
1300
1301 // Sets the mesh dimension.
1302 {
1303 Integer wanted_dimension = mesh_dimension;
1304 wanted_dimension = mesh->parallelMng()->reduce(Parallel::ReduceMax, wanted_dimension);
1305 //mesh->setDimension(wanted_dimension);
1306 mesh_build_info.setMeshDimension(wanted_dimension);
1307 }
1308
1309 mesh_build_info.allocateMesh();
1310
1311 // Positions the coordinates
1312 {
1313 VariableNodeReal3& nodes_coord_var(mesh->nodesCoordinates());
1314 ENUMERATE_NODE (inode, mesh->allNodes()) {
1315 Node node = *inode;
1316 nodes_coord_var[inode] = node_coords[node.uniqueId().asInt32()];
1317 }
1318 }
1319
1320 // Now, check if there is data associated with the file
1321 bool r = _readData(mesh, vtk_file, use_internal_partition, IK_Cell, cells_local_id, nb_node);
1322 debug() << "Reading _readData OK";
1323
1324 return r;
1325}
1326
1327/*---------------------------------------------------------------------------*/
1328/*---------------------------------------------------------------------------*/
1329
1339_readFacesMesh(IMesh* mesh, const String& file_name, const String& dir_name,
1340 bool use_internal_partition)
1341{
1342 ARCANE_UNUSED(dir_name);
1343
1344 std::ifstream ifile(file_name.localstr(), std::ifstream::binary);
1345 if (!ifile) {
1346 info() << "No face descriptor file found '" << file_name << "'";
1347 return;
1348 }
1349
1350 VtkFile vtk_file(&ifile);
1351 const char* buf = 0;
1352
1353 // Reading the description
1354 String title = vtk_file.getNextLine();
1355 info() << "Reading VTK file '" << file_name << "'";
1356 info() << "Title of VTK file: " << title;
1357
1358 String format = vtk_file.getNextLine();
1359 if (VtkFile::isEqualString(format, "BINARY")) {
1360 vtk_file.setIsBinaryFile(true);
1361 }
1362
1363 eMeshType mesh_type = VTK_MT_Unknown;
1364 // Reading the mesh type
1365 // TODO: in parallel, with use_internal_partition true, only processor 0
1366 // reads the data. In this case, it is unnecessary for others to open the file.
1367 {
1368 buf = vtk_file.getNextLine();
1369 std::istringstream mesh_type_line(buf);
1370 std::string dataset_str;
1371 std::string mesh_type_str;
1372 mesh_type_line >> ws >> dataset_str >> ws >> mesh_type_str;
1373 vtk_file.checkString(dataset_str, "DATASET");
1374
1375 if (VtkFile::isEqualString(mesh_type_str, "UNSTRUCTURED_GRID")) {
1376 mesh_type = VTK_MT_UnstructuredGrid;
1377 }
1378
1379 if (mesh_type == VTK_MT_Unknown) {
1380 error() << "Face descriptor file type must be 'UNSTRUCTURED_GRID' (format=" << mesh_type_str << "')";
1381 return;
1382 }
1383 }
1384
1385 {
1386 IParallelMng* pm = mesh->parallelMng();
1387 Integer nb_face = 0;
1388 Int32 sid = pm->commRank();
1389
1390 UniqueArray<Int32> faces_local_id;
1391
1392 // If we use the internal partitioner, only the subdomain reads the mesh
1393 bool need_read = true;
1394 if (use_internal_partition)
1395 need_read = (sid == 0);
1396
1397 if (need_read) {
1398 {
1399 // Reads nodes, but does not keep their coordinates because it is
1400 // not necessary.
1401 UniqueArray<Real3> node_coords;
1402 _readNodesUnstructuredGrid(mesh, vtk_file, node_coords);
1403 //nb_node = node_coords.size();
1404 }
1405
1406 // Reading face info
1407 // Reading connectivity
1408 UniqueArray<Integer> faces_nb_node;
1409 UniqueArray<Int64> faces_connectivity;
1410 UniqueArray<ItemTypeId> faces_type;
1411 _readCellsUnstructuredGrid(mesh, vtk_file, faces_nb_node, faces_type, faces_connectivity);
1412 nb_face = faces_nb_node.size();
1413 //nb_face_node = faces_connectivity.size();
1414
1415 // We must retrieve the localId() of the faces from the connectivity
1416 faces_local_id.resize(nb_face);
1417 {
1418 IMeshUtilities* mu = mesh->utilities();
1419 mu->getFacesLocalIdFromConnectivity(faces_type, faces_connectivity, faces_local_id);
1420 }
1421 }
1422
1423 // Now, check if there is data associated with the files
1424 _readData(mesh, vtk_file, use_internal_partition, IK_Face, faces_local_id, 0);
1425 }
1426}
1427
1428/*---------------------------------------------------------------------------*/
1429/*---------------------------------------------------------------------------*/
1430
1443_readData(IMesh* mesh, VtkFile& vtk_file, bool use_internal_partition,
1444 eItemKind cell_kind, Int32ConstArrayView local_id, Integer nb_node)
1445{
1446 // Only the master subdomain reads the values. However, the other
1447 // subdomains must know the list of variables and groups created.
1448 // If a data item is named 'GROUP_*', it is considered a
1449 // group
1450
1451 IParallelMng* pm = mesh->parallelMng();
1452 IVariableMng* variable_mng = mesh->variableMng();
1453
1454 Int32 sid = pm->commRank();
1455
1456 // If there is no data, return immediately.
1457 {
1458 Byte has_data = 1;
1459 if ((sid == 0) && vtk_file.isEof())
1460 has_data = 0;
1461
1462 ByteArrayView bb(1, &has_data);
1463 pm->broadcast(bb, 0);
1464 // No data.
1465 if (!has_data)
1466 return false;
1467 }
1468
1469 OStringStream created_infos_str;
1470 created_infos_str() << "<?xml version='1.0' ?>\n";
1471 created_infos_str() << "<infos>";
1472
1473 Integer nb_cell_kind = mesh->nbItem(cell_kind);
1474 const char* buf = 0;
1475
1476 if (sid == 0) {
1477 bool reading_node = false;
1478 bool reading_cell = false;
1479 while (((buf = vtk_file.getNextLine()) != 0) && !vtk_file.isEof()) {
1480 debug() << "Read line";
1481 std::istringstream iline(buf);
1482 std::string data_str;
1483 iline >> data_str;
1484
1485 // If we have a "CELL_DATA" block.
1486 if (VtkFile::isEqualString(data_str, "CELL_DATA")) {
1487 Integer nb_item = 0;
1488 iline >> ws >> nb_item;
1489 reading_node = false;
1490 reading_cell = true;
1491 if (nb_item != nb_cell_kind)
1492 error() << "Size expected = " << nb_cell_kind << " found = " << nb_item;
1493 }
1494
1495 // If we have a "POINT_DATA" block.
1496 else if (VtkFile::isEqualString(data_str, "POINT_DATA")) {
1497 Integer nb_item = 0;
1498 iline >> ws >> nb_item;
1499 reading_node = true;
1500 reading_cell = false;
1501 if (nb_item != nb_node)
1502 error() << "Size expected = " << nb_node << " found = " << nb_item;
1503 }
1504
1505 // If we have a "FIELD" block.
1506 else if (VtkFile::isEqualString(data_str, "FIELD")) {
1507 std::string name_str;
1508 int nb_fields;
1509
1510 iline >> ws >> name_str >> ws >> nb_fields;
1511
1512 Integer nb_item = 0;
1513 std::string type_str;
1514 std::string s_name_str;
1515 int nb_component = 1;
1516 bool is_group = false;
1517
1518 for (Integer i = 0; i < nb_fields; i++) {
1519 buf = vtk_file.getNextLine();
1520 std::istringstream iline(buf);
1521 iline >> ws >> s_name_str >> ws >> nb_component >> ws >> nb_item >> ws >> type_str;
1522
1523 if (nb_item != nb_cell_kind && reading_cell && !reading_node)
1524 error() << "Size expected = " << nb_cell_kind << " found = " << nb_item;
1525
1526 if (nb_item != nb_node && !reading_cell && reading_node)
1527 error() << "Size expected = " << nb_node << " found = " << nb_item;
1528
1529 String name_str = s_name_str;
1530 String cstr = name_str.substring(0, 6);
1531
1532 if (cstr == "GROUP_") {
1533 is_group = true;
1534 String new_name = name_str.substring(6);
1535 debug() << "** ** ** GROUP ! name=" << new_name;
1536 name_str = new_name;
1537 }
1538
1539 if (is_group) {
1540 if (!VtkFile::isEqualString(type_str, "int")) {
1541 error() << "Group type must be 'int', found=" << type_str;
1542 return true;
1543 }
1544
1545 if (reading_node) {
1546 created_infos_str() << "<node-group name='" << name_str << "'/>";
1547 _readNodeGroup(mesh, vtk_file, name_str, nb_node);
1548 }
1549
1550 if (reading_cell) {
1551 created_infos_str() << "<cell-group name='" << name_str << "'/>";
1552 _readItemGroup(mesh, vtk_file, name_str, nb_cell_kind, cell_kind, local_id);
1553 }
1554 }
1555
1556 // TODO : See an example if possible.
1557 else {
1558 if (!VtkFile::isEqualString(type_str, "float") && !VtkFile::isEqualString(type_str, "double")) {
1559 error() << "Expecting 'float' or 'double' data type, found=" << type_str;
1560 return true;
1561 }
1562
1563 if (reading_node) {
1564 fatal() << "Unable to read POINT_DATA: feature not implemented";
1565 }
1566 if (reading_cell) {
1567 created_infos_str() << "<cell-variable name='" << name_str << "'/>";
1568
1569 if (cell_kind != IK_Cell)
1570 throw IOException("Unable to read face variables: feature not supported");
1571
1572 _readCellVariable(mesh, vtk_file, name_str, nb_cell_kind);
1573 }
1574 }
1575 }
1576 }
1577
1578 else {
1579 // Reading the values (CELL_DATA or POINT_DATA block)
1580 if (reading_node || reading_cell) {
1581 std::string type_str;
1582 std::string s_name_str;
1583 //String name_str;
1584 bool is_group = false;
1585 int nb_component = 1;
1586
1587 iline >> ws >> s_name_str >> ws >> type_str >> ws >> nb_component;
1588 debug() << "** ** ** READNAME: name=" << s_name_str << " type=" << type_str;
1589
1590 String name_str = s_name_str;
1591 String cstr = name_str.substring(0, 6);
1592
1593 if (cstr == "GROUP_") {
1594 is_group = true;
1595 String new_name = name_str.substring(6);
1596 info() << "** ** ** GROUP ! name=" << new_name;
1597 name_str = new_name;
1598 }
1599
1600 if (!VtkFile::isEqualString(data_str, "SCALARS")) {
1601 error() << "Expecting 'SCALARS' data type, found=" << data_str;
1602 return true;
1603 }
1604
1605 if (is_group) {
1606 if (!VtkFile::isEqualString(type_str, "int")) {
1607 error() << "Group type must be 'int', found=" << type_str;
1608 return true;
1609 }
1610
1611 // To read LOOKUP_TABLE
1612 buf = vtk_file.getNextLine();
1613
1614 if (reading_node) {
1615 created_infos_str() << "<node-group name='" << name_str << "'/>";
1616 _readNodeGroup(mesh, vtk_file, name_str, nb_node);
1617 }
1618
1619 if (reading_cell) {
1620 created_infos_str() << "<cell-group name='" << name_str << "'/>";
1621 _readItemGroup(mesh, vtk_file, name_str, nb_cell_kind, cell_kind, local_id);
1622 }
1623 }
1624 else {
1625 if (!VtkFile::isEqualString(type_str, "float") && !VtkFile::isEqualString(type_str, "double")) {
1626 error() << "Expecting 'float' or 'double' data type, found=" << type_str;
1627 return true;
1628 }
1629
1630 // To read LOOKUP_TABLE
1631 /*buf = */ vtk_file.getNextLine();
1632 if (reading_node) {
1633 fatal() << "Unable to read POINT_DATA: feature not implemented";
1634 }
1635 if (reading_cell) {
1636 created_infos_str() << "<cell-variable name='" << name_str << "'/>";
1637
1638 if (cell_kind != IK_Cell)
1639 throw IOException("Unable to read face variables: feature not supported");
1640
1641 _readCellVariable(mesh, vtk_file, name_str, nb_cell_kind);
1642 }
1643 }
1644 }
1645 else {
1646 error() << "Expecting value CELL_DATA or POINT_DATA, found='" << data_str << "'";
1647 return true;
1648 }
1649 }
1650 }
1651 }
1652 created_infos_str() << "</infos>";
1653 if (use_internal_partition) {
1654 ByteUniqueArray bytes;
1655 if (sid == 0) {
1656 String str = created_infos_str.str();
1657 ByteConstArrayView bv = str.utf8();
1658 Integer len = bv.size();
1659 bytes.resize(len + 1);
1660 bytes.copy(bv);
1661 }
1662
1663 pm->broadcastMemoryBuffer(bytes, 0);
1664
1665 if (sid != 0) {
1666 String str = String::fromUtf8(bytes);
1667 info() << "FOUND STR=" << bytes.size() << " " << str;
1669 XmlNode doc_node = doc->documentNode();
1670
1671 // Reading variables
1672 {
1673 XmlNodeList vars = doc_node.documentElement().children("cell-variable");
1674 for (XmlNode xnode : vars) {
1675 String name = xnode.attrValue("name");
1676 info() << "Building variable: " << name;
1678 variable_mng->_internalApi()->addAutoDestroyVariable(var);
1679 }
1680 }
1681 // Reading cell groups
1682 {
1683 XmlNodeList vars = doc_node.documentElement().children("cell-group");
1684 IItemFamily* cell_family = mesh->itemFamily(cell_kind);
1685 for (XmlNode xnode : vars) {
1686 String name = xnode.attrValue("name");
1687 info() << "Building group: " << name;
1688 cell_family->createGroup(name);
1689 }
1690 }
1691
1692 // Reading node groups
1693 {
1694 XmlNodeList vars = doc_node.documentElement().children("node-group");
1695 IItemFamily* node_family = mesh->nodeFamily();
1696 for (XmlNode xnode : vars) {
1697 String name = xnode.attrValue("name");
1698 info() << "Create node group: " << name;
1699 node_family->createGroup(name);
1700 }
1701 }
1702 }
1703 }
1704 return false;
1705}
1706
1707/*---------------------------------------------------------------------------*/
1708/*---------------------------------------------------------------------------*/
1709
1719_createFaceGroup(IMesh* mesh, const String& name, Int32ConstArrayView faces_lid)
1720{
1721 info() << "Building face group '" << name << "'"
1722 << " size=" << faces_lid.size();
1723
1724 mesh->faceFamily()->createGroup(name, faces_lid);
1725}
1726
1727/*---------------------------------------------------------------------------*/
1728/*---------------------------------------------------------------------------*/
1729
1739_readCellVariable(IMesh* mesh, VtkFile& vtk_file, const String& var_name, Integer nb_cell)
1740{
1741 //TODO Perform the correct conversion from uniqueId() to localId()
1742 info() << "Reading values for variable: " << var_name << " n=" << nb_cell;
1743 auto* var = new VariableCellReal(VariableBuildInfo(mesh, var_name));
1744 mesh->variableMng()->_internalApi()->addAutoDestroyVariable(var);
1745 RealArrayView values(var->asArray());
1746 for (Integer i = 0; i < nb_cell; ++i) {
1747 Real v = vtk_file.getDouble();
1748 values[i] = v;
1749 }
1750 _readMetadata(mesh, vtk_file);
1751 info() << "Variable build finished: " << vtk_file.isEof();
1752}
1753
1754/*---------------------------------------------------------------------------*/
1755/*---------------------------------------------------------------------------*/
1756
1768_readItemGroup(IMesh* mesh, VtkFile& vtk_file, const String& name, Integer nb_item,
1769 eItemKind ik, Int32ConstArrayView local_id)
1770{
1771 IItemFamily* item_family = mesh->itemFamily(ik);
1772 info() << "Reading group info for group: " << name;
1773
1774 Int32UniqueArray ids;
1775 for (Integer i = 0; i < nb_item; ++i) {
1776 Integer v = vtk_file.getInt();
1777 if (v != 0)
1778 ids.add(local_id[i]);
1779 }
1780 info() << "Building group: " << name << " nb_element=" << ids.size();
1781
1782 item_family->createGroup(name, ids);
1783
1784 _readMetadata(mesh, vtk_file);
1785}
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() << "Reading node group info for group: " << 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() << "Creating group: " << 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/*---------------------------------------------------------------------------*/
1856
1870writeMeshToFile(IMesh* mesh, const String& file_name)
1871{
1872 String fname = file_name;
1873 // Append the '.vtk' extension if it is not present.
1874 if (!fname.endsWith(".vtk"))
1875 fname = fname + ".vtk";
1876 _writeMeshToFile(mesh, fname, IK_Cell);
1877 _writeMeshToFile(mesh, fname + "faces.vtk", IK_Face);
1878 return false;
1879}
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 // Save nodes
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 // Save cells or 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 // The type must be consistent with vtkCellType.h
1941 ofile << "CELL_TYPES " << nb_cell_kind << "\n";
1942 ENUMERATE_ (ItemWithNodes, iitem, cell_kind_family->allItems()) {
1943 const ItemTypeInfo* iti = iitem->typeInfo();
1944 // Check if the type is a polygon for VTK (dimension 2 and more than 4 nodes)
1945 int type = VTK_BAD_ARCANE_TYPE;
1946 if (iti->isPolygon())
1947 type = VTK_POLYGON;
1948 else
1949 type = arcaneToVtkCellType(iti);
1950 ofile << type << '\n';
1951 }
1952 }
1953
1954 // If we are in the cell mesh, save node groups.
1955 if (cell_kind == IK_Cell) {
1956 ofile << "POINT_DATA " << nb_node << "\n";
1957 _saveGroups(mesh->itemFamily(IK_Node), ofile);
1958 }
1959
1960 // Save cell groups
1961 ofile << "CELL_DATA " << nb_cell_kind << "\n";
1962 _saveGroups(mesh->itemFamily(cell_kind), ofile);
1963}
1964
1965/*---------------------------------------------------------------------------*/
1966/*---------------------------------------------------------------------------*/
1967
1968void VtkLegacyMeshWriter::
1969_saveGroups(IItemFamily* family, std::ostream& ofile)
1970{
1971 info() << "Saving groups for family name=" << family->name();
1972 UniqueArray<char> in_group_list(family->maxLocalId());
1973 for (ItemGroupCollection::Enumerator igroup(family->groups()); ++igroup;) {
1974 ItemGroup group = *igroup;
1975 // No need to save the group of all entities
1976 if (group == family->allItems())
1977 continue;
1978 //HACK: to be removed
1979 if (group.name() == "OuterFaces")
1980 continue;
1981 ofile << "SCALARS GROUP_" << group.name() << " int 1\n";
1982 ofile << "LOOKUP_TABLE default\n";
1983 in_group_list.fill('0');
1984 ENUMERATE_ITEM (iitem, group) {
1985 in_group_list[(*iitem).localId()] = '1';
1986 }
1987 ENUMERATE_ITEM (iitem, family->allItems()) {
1988 ofile << in_group_list[(*iitem).localId()] << '\n';
1989 }
1990 }
1991}
1992
1993/*---------------------------------------------------------------------------*/
1994/*---------------------------------------------------------------------------*/
1995
1996/*---------------------------------------------------------------------------*/
1997/*---------------------------------------------------------------------------*/
1998
1999class VtkMeshReader
2000: public AbstractService
2001, public IMeshReader
2002{
2003 public:
2004
2005 explicit VtkMeshReader(const ServiceBuildInfo& sbi)
2006 : AbstractService(sbi)
2007 {}
2008
2009 public:
2010
2011 bool allowExtension(const String& str) override { return str == "vtk"; }
2012 eReturnType readMeshFromFile(IPrimaryMesh* mesh, const XmlNode& mesh_node, const String& file_name,
2013 const String& dir_name, bool use_internal_partition) override
2014
2015 {
2016 ARCANE_UNUSED(mesh_node);
2017 VtkMeshIOService vtk_service(traceMng());
2018 bool ret = vtk_service.readMesh(mesh, file_name, dir_name, use_internal_partition);
2019 if (ret)
2020 return RTError;
2021
2022 return RTOk;
2023 }
2024};
2025
2026/*---------------------------------------------------------------------------*/
2027/*---------------------------------------------------------------------------*/
2028
2030
2032 ServiceProperty("VtkLegacyMeshReader", ST_SubDomain),
2034
2035/*---------------------------------------------------------------------------*/
2036/*---------------------------------------------------------------------------*/
2037
2038/*---------------------------------------------------------------------------*/
2039/*---------------------------------------------------------------------------*/
2040
2041class VtkLegacyCaseMeshReader
2042: public AbstractService
2043, public ICaseMeshReader
2044{
2045 public:
2046
2047 class Builder
2048 : public IMeshBuilder
2049 {
2050 public:
2051
2052 explicit Builder(ITraceMng* tm, const CaseMeshReaderReadInfo& read_info)
2053 : m_trace_mng(tm)
2054 , m_read_info(read_info)
2055 {}
2056
2057 public:
2058
2059 void fillMeshBuildInfo(MeshBuildInfo& build_info) override
2060 {
2061 ARCANE_UNUSED(build_info);
2062 }
2064 {
2065 VtkMeshIOService vtk_service(m_trace_mng);
2066 String fname = m_read_info.fileName();
2067 m_trace_mng->info() << "VtkLegacy Reader (ICaseMeshReader) file_name=" << fname;
2068 bool ret = vtk_service.readMesh(pm, fname, m_read_info.directoryName(), m_read_info.isParallelRead());
2069 if (ret)
2070 ARCANE_FATAL("Can not read VTK File");
2071 }
2072
2073 private:
2074
2075 ITraceMng* m_trace_mng;
2076 CaseMeshReaderReadInfo m_read_info;
2077 };
2078
2079 public:
2080
2081 explicit VtkLegacyCaseMeshReader(const ServiceBuildInfo& sbi)
2082 : AbstractService(sbi)
2083 {}
2084
2085 public:
2086
2088 {
2089 IMeshBuilder* builder = nullptr;
2090 if (read_info.format() == "vtk")
2091 builder = new Builder(traceMng(), read_info);
2092 return makeRef(builder);
2093 }
2094};
2095
2096/*---------------------------------------------------------------------------*/
2097/*---------------------------------------------------------------------------*/
2098
2100 ServiceProperty("VtkLegacyCaseMeshReader", ST_SubDomain),
2102
2103/*---------------------------------------------------------------------------*/
2104/*---------------------------------------------------------------------------*/
2105
2106} // End namespace Arcane
2107
2108/*---------------------------------------------------------------------------*/
2109/*---------------------------------------------------------------------------*/
#define ARCANE_THROW(exception_class,...)
Macro for throwing an exception with formatting.
#define ARCANE_FATAL(...)
Macro throwing a FatalErrorException.
Types and macros for iterating over mesh entities.
#define ENUMERATE_FACE(name, group)
Generic enumerator for a face group.
#define ENUMERATE_(type, name, group)
Generic enumerator for an entity group.
#define ENUMERATE_ITEM(name, group)
Generic enumerator for a node group.
#define ENUMERATE_NODE(name, group)
Generic enumerator for a node group.
#define ARCANE_REGISTER_SUB_DOMAIN_FACTORY(aclass, ainterface, aname)
Registers a factory service for the class aclass.
#define ARCANE_SERVICE_INTERFACE(ainterface)
Macro to declare an interface when registering a service.
Integer size() const
Number of elements in the vector.
Base class of a service.
AbstractService(const ServiceBuildInfo &)
Constructor from a ServiceBuildInfo.
Base class for 1D data vectors.
void fill(ConstReferenceType value)
Fills the array with the value value.
void resize(Int64 s)
Changes the number of elements in the array to s.
void add(ConstReferenceType val)
Adds element val to the end of the array.
ArrayView< T > subView(Int64 abegin, Integer asize)
Sub-view starting from element abegin and containing asize elements.
Necessary information for reading a mesh file.
Constant view of an array of type T.
constexpr Integer size() const noexcept
Number of elements in the array.
Face of a cell.
Definition Item.h:1032
Information about the floating-point type.
Definition Limits.h:49
Interface for the mesh reading service from the dataset.
Interface of an entity family.
Definition IItemFamily.h:83
virtual ItemGroupCollection groups() const =0
Collection of groups in this family.
virtual ItemGroup allItems() const =0
Group of all entities.
virtual ItemGroup createGroup(const String &name, Int32ConstArrayView local_ids, bool do_override=false)=0
Creates an entity group named name containing the entities local_ids.
virtual Int32 maxLocalId() const =0
virtual String name() const =0
Family name.
virtual Integer nbItem() const =0
Number of entities.
Interface of a mesh creation/reading service.
Interface of the service managing the reading of a mesh.
Definition IMeshReader.h:33
eReturnType
Types of return codes for a read or write operation.
Definition IMeshReader.h:38
@ RTError
Error during the operation.
Definition IMeshReader.h:40
@ RTOk
Operation successfully performed.
Definition IMeshReader.h:39
Interface of a class providing utility functions on meshes.
virtual void getFacesLocalIdFromConnectivity(ConstArrayView< ItemTypeId > items_type, ConstArrayView< Int64 > items_connectivity, ArrayView< Int32 > local_ids, bool allow_null=false)=0
Searches for the local IDs of faces based on their connectivity.
Interface of a mesh writing service.
Definition IMeshWriter.h:38
Exception when an input/output error is detected.
Definition IOException.h:34
Interface of the parallelism manager for a subdomain.
virtual Int32 commRank() const =0
Rank of this instance in the communicator.
virtual void broadcastMemoryBuffer(ByteArray &bytes, Int32 rank)=0
Performs a broadcast of a memory region.
virtual void addAutoDestroyVariable(VariableRef *var)=0
Adds the variable to the list of variables that are kept until the end of execution.
Variable manager interface.
virtual IVariableMngInternal * _internalApi()=0
Internal Arcane API.
static IXmlDocumentHolder * loadFromBuffer(Span< const Byte > buffer, const String &name, ITraceMng *tm)
Loads an XML document.
Definition DomUtils.cc:425
Mesh entity group.
Definition ItemGroup.h:51
const String & name() const
Group name.
Definition ItemGroup.h:81
Type of an entity (Item).
Definition ItemTypeId.h:33
Info on a mesh entity type.
Int16 dimension() const
Dimension of the element (<0 if unknown).
bool isPolygon() const
Indicates if the type is a polygon.
Mesh entity type manager.
Definition ItemTypeMng.h:66
ItemTypeInfo * typeFromId(Integer id) const
Type corresponding to the number id.
Mesh element based on nodes (Edge,Face,Cell).
Definition Item.h:773
NodeConnectedListViewType nodes() const
List of nodes of the entity.
Definition Item.h:843
Int32 nbNode() const
Number of nodes of the entity.
Definition Item.h:837
constexpr Int32 localId() const
Local identifier of the entity in the processor subdomain.
Definition Item.h:233
ItemUniqueId uniqueId() const
Unique identifier across all domains.
Definition Item.h:239
Parameters necessary for building a mesh.
Node of a mesh.
Definition Item.h:598
Output stream linked to a String.
Class managing a 3-dimensional real vector.
Definition Real3.h:132
Reference to an instance.
Encapsulation of an automatically destructing pointer.
Definition ScopedPtr.h:44
Structure containing the information to create a service.
Service creation properties.
String lower() const
Transforms all characters in the string to lowercase.
Definition String.cc:480
const char * localstr() const
Returns the conversion of the instance into UTF-8 encoding.
Definition String.cc:229
ByteConstArrayView utf8() const
Returns the conversion of the instance into UTF-8 encoding.
Definition String.cc:277
bool endsWith(const String &s) const
Indicates if the string ends with the characters of s.
Definition String.cc:1095
String substring(Int64 pos) const
Substring starting at position pos.
Definition String.cc:1126
TraceAccessor(ITraceMng *m)
Constructs an accessor via the trace manager m.
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flow for a debug message.
TraceMessage fatal() const
Flow for a fatal error message.
TraceMessage info() const
Flow for an information message.
TraceMessage error() const
Flow for an error message.
ITraceMng * traceMng() const
Trace manager.
1D data vector with value semantics (STL style).
Information for allocating entities of an unstructured mesh.
void allocateMesh()
Allocates the mesh with the cells added during the call to addCell().
void preAllocate(Int32 nb_cell, Int64 nb_connectivity_node)
Pre-allocate the memory.
void addCell(ItemTypeId type_id, Int64 cell_uid, SmallSpan< const Int64 > nodes_uid)
Adds a cell to the mesh.
Parameters necessary for building a variable.
const char * getCurrentLine()
Allows returning the line present in the buffer.
bool isEmptyNextLine()
Allows checking if the next line is empty.
bool m_is_init
Has at least one line been read.
std::istream * m_stream
The stream.
static bool isEqualString(const String &current_value, const String &expected_value)
Allows checking if expected_value == current_value.
double getDouble()
Allows retrieving the next double.
bool m_is_eof
Is the end of the file reached.
char m_buf[BUFSIZE]
The buffer containing the read line.
void getBinary(T &type)
Allows retrieving the next binary number.
int getInt()
Allows retrieving the next integer.
void checkString(const String &current_value, const String &expected_value)
Allows checking if expected_value == current_value.
bool m_is_binary_file
Is this a file containing binary data.
bool m_need_reread_current_line
Should the same line be reread.
const char * getNextLine()
Allows retrieving the next line from the file.
float getFloat()
Allows retrieving the following float.
void allocateMeshItems(IPrimaryMesh *pm) override
Allocates the mesh entities managed by this service.
void fillMeshBuildInfo(MeshBuildInfo &build_info) override
Fills build_info with the necessary information to create the mesh.
Ref< IMeshBuilder > createBuilder(const CaseMeshReaderReadInfo &read_info) const override
Returns a builder to create and read the mesh whose information is specified in read_info.
bool writeMeshToFile(IMesh *mesh, const String &file_name) override
Writing the mesh to VTK.
void _writeMeshToFile(IMesh *mesh, const String &file_name, eItemKind cell_kind)
Writes the mesh in VTK format.
void build() override
Build-level construction of the service.
Mesh file reader for legacy VTK format.
void _createFaceGroup(IMesh *mesh, const String &name, Int32ConstArrayView faces_lid)
Allows creating a face group named "name" composed of faces having the IDs included in "faces_lid".
bool _readMetadata(IMesh *mesh, VtkFile &vtk_file)
Read metadata.
void _readCellVariable(IMesh *mesh, VtkFile &vtk_file, const String &name_str, Integer nb_cell)
Allows creating a cell variable from the information in the vtk file.
bool _readStructuredGrid(IPrimaryMesh *mesh, VtkFile &, bool use_internal_partition)
Allows reading a vtk file containing a STRUCTURED_GRID.
void _readNodesUnstructuredGrid(IMesh *mesh, VtkFile &vtk_file, Array< Real3 > &node_coords)
Read nodes and their coordinates.
void _readItemGroup(IMesh *mesh, VtkFile &vtk_file, const String &name_str, Integer nb_item, eItemKind ik, ConstArrayView< Int32 > local_id)
Allows creating an item group.
void _readCellsUnstructuredGrid(IMesh *mesh, VtkFile &vtk_file, Array< Int32 > &cells_nb_node, Array< ItemTypeId > &cells_type, Array< Int64 > &cells_connectivity)
Read cells and their connectivity.
void _readFacesMesh(IMesh *mesh, const String &file_name, const String &dir_name, bool use_internal_partition)
Allows reading the truc.vtkfaces.vtk file (if it exists).
bool _readData(IMesh *mesh, VtkFile &vtk_file, bool use_internal_partition, eItemKind cell_kind, Int32ConstArrayView local_id, Integer nb_node)
Allows reading supplementary data (POINT_DATA / CELL_DATA).
bool readMesh(IPrimaryMesh *mesh, const String &file_name, const String &dir_name, bool use_internal_partition)
Allows starting the reading of a vtk file.
void _readNodeGroup(IMesh *mesh, VtkFile &vtk_file, const String &name, Integer nb_item)
Allows creating a node group.
bool _readUnstructuredGrid(IPrimaryMesh *mesh, VtkFile &vtk_file, bool use_internal_partition)
Allows reading a vtk file containing an UNSTRUCTURED_GRID.
eReturnType readMeshFromFile(IPrimaryMesh *mesh, const XmlNode &mesh_node, const String &file_name, const String &dir_name, bool use_internal_partition) override
Reads a mesh from a file.
bool allowExtension(const String &str) override
Checks if the service supports files with the extension str.
List of nodes of a DOM tree.
Definition XmlNodeList.h:36
Node of a DOM tree.
Definition XmlNode.h:51
XmlNode documentElement() const
Returns the document element.
Definition XmlNode.cc:565
XmlNodeList children(const String &name) const
Set of child nodes of this node having the name name.
Definition XmlNode.cc:102
#define ARCANE_REGISTER_SERVICE(aclass, a_service_property,...)
Macro for registering a service.
MeshVariableScalarRefT< Cell, Real > VariableCellReal
Real type quantity at cell center.
MeshVariableScalarRefT< Node, Real3 > VariableNodeReal3
Coordinate type quantity at node.
@ ReduceMax
Maximum of values.
-- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature --
ArrayView< Byte > ByteArrayView
C equivalent of a 1D array of characters.
Definition UtilsTypes.h:447
std::int64_t Int64
Signed integer type of 64 bits.
Int32 Integer
Type representing an integer.
ConstArrayView< Int32 > Int32ConstArrayView
C equivalent of a 1D array of 32-bit integers.
Definition UtilsTypes.h:482
@ ST_SubDomain
The service is used at the subdomain level.
UniqueArray< Byte > ByteUniqueArray
Dynamic 1D array of characters.
Definition UtilsTypes.h:335
UniqueArray< Int32 > Int32UniqueArray
Dynamic 1D array of 32-bit integers.
Definition UtilsTypes.h:341
eItemKind
Mesh entity type.
@ IK_Node
Node mesh entity.
@ IK_Cell
Cell mesh entity.
@ IK_Face
Face mesh entity.
std::int16_t Int16
Signed integer type of 16 bits.
double Real
Type representing a real number.
ConstArrayView< Byte > ByteConstArrayView
C equivalent of a 1D array of characters.
Definition UtilsTypes.h:476
unsigned char Byte
Type of a byte.
Definition BaseTypes.h:43
auto makeRef(InstanceType *t) -> Ref< InstanceType >
Creates a reference on a pointer.
ArrayView< Real > RealArrayView
C equivalent of a 1D array of reals.
Definition UtilsTypes.h:459
std::int32_t Int32
Signed integer type of 32 bits.
Real y
second component of the triplet
Definition Real3.h:36
Real z
third component of the triplet
Definition Real3.h:37
Real x
first component of the triplet
Definition Real3.h:35