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"
27#include "arcane/core/BasicService.h"
28#include "arcane/core/FactoryService.h"
29#include "arcane/core/ICaseMeshReader.h"
30#include "arcane/core/IItemFamily.h"
31#include "arcane/core/IPrimaryMesh.h"
32#include "arcane/core/IMeshBuilder.h"
33#include "arcane/core/IMeshReader.h"
34#include "arcane/core/IMeshUtilities.h"
35#include "arcane/core/IMeshWriter.h"
36#include "arcane/core/IParallelMng.h"
37#include "arcane/core/IVariableAccessor.h"
38#include "arcane/core/IXmlDocumentHolder.h"
39#include "arcane/core/Item.h"
41#include "arcane/core/IVariableMng.h"
42#include "arcane/core/VariableTypes.h"
43#include "arcane/core/XmlNode.h"
44#include "arcane/core/XmlNodeList.h"
45#include "arcane/core/UnstructuredMeshAllocateBuildInfo.h"
47#include "arcane/core/internal/IVariableMngInternal.h"
49#include "arcane/std/internal/VtkCellTypes.h"
61using namespace Arcane::VtkUtils;
120 VTK_MT_StructuredGrid,
121 VTK_MT_UnstructuredGrid
173 static const int BUFSIZE = 10000;
177 explicit VtkFile(std::istream* stream)
271 throw IOException(
"VtkFile::isEmptyNextLine()",
"Unexpected EndOfFile");
300 for (
int i = 0; i < BUFSIZE &&
m_buf[i] !=
'\0'; ++i) {
304 if (
m_buf[i] ==
'#') {
313 for (
int i = 0; i < BUFSIZE &&
m_buf[i] !=
'\0'; ++i) {
314 if (
m_buf[i] ==
'\r') {
330 throw IOException(
"VtkFile::isEmptyNextLine()",
"Not Good");
355 throw IOException(
"VtkFile::isEmptyNextLine()",
"Unexpected EndOfFile");
378 for (
int i = 0; i < BUFSIZE &&
m_buf[i] !=
'\0'; ++i) {
382 if (
m_buf[i] ==
'#') {
391 for (
int i = 0; i < BUFSIZE &&
m_buf[i] !=
'\0'; ++i) {
392 if (
m_buf[i] ==
'\r') {
400 throw IOException(
"VtkFile::getNextLine()",
"Not good");
419 (*m_stream) >> ws >> v;
424 throw IOException(
"VtkFile::getFloat()",
"Bad float");
443 (*m_stream) >> ws >> v;
448 throw IOException(
"VtkFile::getDouble()",
"Bad double");
467 (*m_stream) >> ws >> v;
487 constexpr size_t sizeofT =
sizeof(T);
497 for (
size_t i = 0; i <
sizeofT; i++) {
610 info() <<
"Titre du fichier VTK : " << title.
localstr();
646 error() <<
"Support exists only for 'STRUCTURED_GRID' and 'UNSTRUCTURED_GRID' formats (format=" <<
mesh_type_str <<
"')";
650 debug() <<
"Lecture en-tête OK";
654 case VTK_MT_StructuredGrid:
658 case VTK_MT_UnstructuredGrid:
660 debug() <<
"Lecture _readUnstructuredGrid OK";
664 debug() <<
"Lecture _readFacesMesh OK";
694 const char* buf =
nullptr;
700 std::istringstream
iline(buf);
705 error() <<
"Syntax error while reading grid dimensions";
722 std::istringstream
iline(buf);
727 error() <<
"Syntax error while reading grid dimensions";
767 for (Integer x = 0; x <
nb_node_x; ++x) {
768 for (Integer z = 0; z <
nb_node_z; ++z) {
769 for (Integer y = 0; y <
nb_node_y; ++y) {
801 for (Integer z = 0; z <
nb_cell_z; ++z) {
802 for (Integer y = 0; y <
nb_cell_y; ++y) {
803 for (Integer x = 0; x <
nb_cell_x; ++x) {
842 for (Integer z = 0; z <
nb_node_z; ++z) {
843 for (Integer y = 0; y <
nb_node_y; ++y) {
844 for (Integer x = 0; x <
nb_node_x; ++x) {
855 for (Integer z = 0; z <
nb_node_z; ++z) {
856 for (Integer y = 0; y <
nb_node_y; ++y) {
857 for (Integer x = 0; x <
nb_node_x; ++x) {
868 for (Integer z = 0; z <
nb_node_z; ++z) {
869 for (Integer y = 0; y <
nb_node_y; ++y) {
870 for (Integer x = 0; x <
nb_node_x; ++x) {
881 throw IOException(
"_readStructuredGrid",
"Invalid type name");
903 Integer face_local_id = face.
localId();
974 const char* func_name =
"VtkMeshIOService::_readNodesUnstructuredGrid()";
975 const char* buf =
vtk_file.getNextLine();
976 std::istringstream
iline(buf);
984 throw IOException(func_name,
"Syntax error while reading number of nodes");
989 throw IOException(A_FUNCINFO, String::format(
"Invalid number of nodes: n={0}",
nb_node));
997 for (Integer i = 0; i <
nb_node; ++i) {
1005 for (Integer i = 0; i <
nb_node; ++i) {
1013 for (Integer i = 0; i <
nb_node; ++i) {
1021 throw IOException(func_name,
"Invalid type name");
1047 ARCANE_UNUSED(mesh);
1049 const char* func_name =
"VtkMeshIOService::_readCellsUnstructuredGrid()";
1050 const char* buf =
vtk_file.getNextLine();
1055 std::istringstream
iline(buf);
1058 Integer nb_cell_node = 0;
1063 throw IOException(func_name,
"Syntax error while reading cells");
1067 if (
nb_cell < 0 || nb_cell_node < 0) {
1069 String::format(
"Invalid dimensions: nb_cell={0} nb_cell_node={1}",
1073 cells_nb_node.resize(
nb_cell);
1075 cells_connectivity.resize(nb_cell_node);
1079 for (Integer i = 0; i <
nb_cell; ++i) {
1081 cells_nb_node[i] = n;
1082 for (Integer
j = 0;
j < n; ++
j) {
1095 std::istringstream
iline(buf);
1101 throw IOException(func_name,
"Syntax error while reading cell types");
1111 for (Integer i = 0; i <
nb_cell; ++i) {
1113 Int16
it = vtkToArcaneCellType(
vtk_ct, cells_nb_node[i]);
1132 ARCANE_UNUSED(mesh);
1242 Integer nb_cell_node = 0;
1261 debug() <<
"Lecture _readNodesUnstructuredGrid OK";
1270 debug() <<
"Lecture _readCellsUnstructuredGrid OK";
1272 nb_cell = cells_nb_node.size();
1273 nb_cell_node = cells_connectivity.size();
1281 for (Integer i = 0; i <
nb_cell; ++i) {
1287 Int16
cell_dim =
itm->typeFromId(cells_type[i])->dimension();
1298 for (Int32 i = 0; i < 4; ++i)
1304 ARCANE_FATAL(
"The mesh contains cells of different dimension. nb0={0} nb1={1} nb2={2} nb3={3}",
1329 debug() <<
"Lecture _readData OK";
1351 std::ifstream
ifile(
file_name.localstr(), std::ifstream::binary);
1353 info() <<
"No face descriptor file found '" <<
file_name <<
"'";
1358 const char* buf = 0;
1362 info() <<
"Titre du fichier VTK : " << title;
1386 error() <<
"Face descriptor file type must be 'UNSTRUCTURED_GRID' (format=" <<
mesh_type_str <<
"')";
1469 pm->broadcast(
bb, 0);
1480 const char* buf = 0;
1486 debug() <<
"Read line";
1487 std::istringstream
iline(buf);
1524 for (Integer i = 0; i <
nb_fields; i++) {
1526 std::istringstream
iline(buf);
1538 if (
cstr ==
"GROUP_") {
1565 error() <<
"Expecting 'float' or 'double' data type, found=" <<
type_str;
1570 fatal() <<
"Unable to read POINT_DATA: feature not implemented";
1576 throw IOException(
"Unable to read face variables: feature not supported");
1599 if (
cstr ==
"GROUP_") {
1607 error() <<
"Expecting 'SCALARS' data type, found=" <<
data_str;
1632 error() <<
"Expecting 'float' or 'double' data type, found=" <<
type_str;
1639 fatal() <<
"Unable to read POINT_DATA: feature not implemented";
1645 throw IOException(
"Unable to read face variables: feature not supported");
1652 error() <<
"Expecting value CELL_DATA or POINT_DATA, found='" <<
data_str <<
"'";
1664 Integer len =
bv.size();
1672 String str = String::fromUtf8(bytes);
1673 info() <<
"FOUND STR=" << bytes.
size() <<
" " << str;
1682 info() <<
"Building variable: " << name;
1694 info() <<
"Building group: " << name;
1705 info() <<
"Create node group: " << name;
1727 info() <<
"Building face group '" << name <<
"'"
1750 mesh->
variableMng()->_internalApi()->addAutoDestroyVariable(
var);
1752 for (Integer i = 0; i <
nb_cell; ++i) {
1757 info() <<
"Variable build finished: " <<
vtk_file.isEof();
1778 info() <<
"Reading group info for group: " << name;
1781 for (Integer i = 0; i <
nb_item; ++i) {
1784 ids.add(local_id[i]);
1786 info() <<
"Building group: " << name <<
" nb_element=" << ids.size();
1808 info() <<
"Lecture infos groupes de noeuds pour le groupe: " << name;
1811 for (Integer i = 0; i <
nb_item; ++i) {
1816 info() <<
"Création groupe: " << name <<
" nb_element=" << ids.size();
1897 throw IOException(
"VtkMeshIOService::writeMeshToFile(): Unable to open file");
1898 ofile <<
"# vtk DataFile Version 2.0\n";
1899 ofile <<
"Maillage Arcane\n";
1901 ofile <<
"DATASET UNSTRUCTURED_GRID\n";
1934 Integer item_nb_node = item.
nbNode();
1935 ofile << item_nb_node;
1946 ofile << type <<
'\n';
1964void VtkLegacyMeshWriter::
1967 info() <<
"Saving groups for family name=" << family->
name();
1975 if (group.
name() ==
"OuterFaces")
1977 ofile <<
"SCALARS GROUP_" << group.
name() <<
" int 1\n";
1978 ofile <<
"LOOKUP_TABLE default\n";
2063 m_trace_mng->
info() <<
"VtkLegacy Reader (ICaseMeshReader) file_name=" <<
fname;
2064 bool ret =
vtk_service.readMesh(pm,
fname, m_read_info.directoryName(), m_read_info.isParallelRead());
2088 return makeRef(builder);
#define ARCANE_THROW(exception_class,...)
Macro pour envoyer une exception avec formattage.
#define ARCANE_FATAL(...)
Macro envoyant une exception FatalErrorException.
#define ARCANE_REGISTER_SUB_DOMAIN_FACTORY(aclass, ainterface, aname)
Enregistre un service de fabrique pour la classe aclass.
#define ARCANE_SERVICE_INTERFACE(ainterface)
Macro pour déclarer une interface lors de l'enregistrement d'un service.
Classe de base d'un service.
Classe de base de service lié à un sous-domaine.
Informations nécessaires pour la lecture d'un fichier de maillage.
Informations sur le type flottant.
Interface du service de lecture du maillage à partir du jeu de données.
Interface d'une famille d'entités.
virtual ItemGroupCollection groups() const =0
Liste des groupes de cette famille.
virtual ItemGroup allItems() const =0
Groupe de toutes les entités.
virtual Int32 maxLocalId() const =0
virtual String name() const =0
Nom de la famille.
virtual IItemFamily * nodeFamily()=0
Retourne la famille des noeuds.
virtual FaceGroup allFaces()=0
Groupe de toutes les faces.
virtual IItemFamily * itemFamily(eItemKind ik)=0
Retourne la famille d'entité de type ik.
virtual Integer nbNode()=0
Nombre de noeuds du maillage.
virtual Integer nbItem(eItemKind ik)=0
Nombre d'éléments du genre ik.
virtual IItemFamily * faceFamily()=0
Retourne la famille des faces.
virtual NodeGroup allNodes()=0
Groupe de tous les noeuds.
Interface d'un service de création/lecture du maillage.
Interface du service gérant la lecture d'un maillage.
eReturnType
Types des codes de retour d'une lecture ou écriture.
@ RTError
Erreur lors de l'opération.
@ RTOk
Opération effectuée avec succès.
Interface d'une classe proposant des fonctions utilitaires sur maillage.
Interface d'un service d'écriture d'un maillage.
virtual IParallelMng * parallelMng()=0
Gestionnaire de parallèlisme.
virtual IMeshUtilities * utilities()=0
Interface des fonctions utilitaires associée.
virtual ItemTypeMng * itemTypeMng() const =0
Gestionnaire de types d'entités associé
virtual IPrimaryMesh * toPrimaryMesh()=0
Retourne l'instance sous la forme d'un IPrimaryMesh.
virtual IVariableMng * variableMng() const =0
Gestionnaire de variable associé
Exception lorsqu'une erreur d'entrée/sortie est détectée.
Interface du gestionnaire de parallélisme pour un sous-domaine.
virtual Int32 commRank() const =0
Rang de cette instance dans le communicateur.
virtual void broadcastMemoryBuffer(ByteArray &bytes, Int32 rank)=0
Effectue un broadcast d'une zone mémoire.
virtual VariableNodeReal3 & nodesCoordinates()=0
Coordonnées des noeuds.
virtual void allocateCells(Integer nb_cell, Int64ConstArrayView cells_infos, bool one_alloc=true)=0
Allocation d'un maillage.
virtual void endAllocate()=0
Indique une fin d'allocation de mailles.
virtual void setDimension(Integer dim)=0
Positionne la dimension du maillage (1D, 2D ou 3D).
Interface du gestionnaire de variables.
static IXmlDocumentHolder * loadFromBuffer(Span< const Byte > buffer, const String &name, ITraceMng *tm)
Charge un document XML.
Groupe d'entités de maillage.
const String & name() const
Nom du groupe.
Type d'une entité (Item).
Gestionnaire des types d'entités de maillage.
Elément de maillage s'appuyant sur des noeuds (Edge,Face,Cell).
NodeConnectedListViewType nodes() const
Liste des noeuds de l'entité
Int32 nbNode() const
Nombre de noeuds de l'entité
constexpr Int32 localId() const
Identifiant local de l'entité dans le sous-domaine du processeur.
ItemUniqueId uniqueId() const
Identifiant unique sur tous les domaines.
Lecteur des fichiers de maillage via la bibliothèque LIMA.
Paramètres nécessaires à la construction d'un maillage.
Flot de sortie lié à une String.
Classe gérant un vecteur de réel de dimension 3.
Structure contenant les informations pour créer un service.
Propriétés de création d'un service.
Informations pour allouer les entités d'un maillage non structuré.
Paramètres nécessaires à la construction d'une variable.
const char * getCurrentLine()
Permet de retourner la ligne présente dans le buffer.
bool isEmptyNextLine()
Permet de voir si la prochaine ligne est vide.
bool m_is_init
Y'a-t-il eu au moins une ligne lue.
std::istream * m_stream
Le stream.
static bool isEqualString(const String ¤t_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 ¤t_value, const String &expected_value)
Permet de vérifier si expected_value == current_value.
bool m_is_binary_file
Est-ce un fichier contenant des données en binaire.
bool m_need_reread_current_line
Doit-on relire la même ligne.
const char * getNextLine()
Permet de récupérer la prochaine ligne du fichier.
float getFloat()
Permet de récupérer le float qui suit.
void allocateMeshItems(IPrimaryMesh *pm) override
Alloue les entités du maillage géré par ce service.
void fillMeshBuildInfo(MeshBuildInfo &build_info) override
Remplit build_info avec les informations nécessaires pour créer le maillage.
Ref< IMeshBuilder > createBuilder(const CaseMeshReaderReadInfo &read_info) const override
Retourne un builder pour créer et lire le maillage dont les informations sont spécifiées dans read_in...
bool writeMeshToFile(IMesh *mesh, const String &file_name) override
Ecriture du maillage au Vtk.
void _writeMeshToFile(IMesh *mesh, const String &file_name, eItemKind cell_kind)
Ecrit le maillage au format Vtk.
void build() override
Construction de niveau build du service.
Lecteur des fichiers de maillage au format Vtk historique (legacy).
void _createFaceGroup(IMesh *mesh, const String &name, Int32ConstArrayView faces_lid)
Permet de créer un groupe de face de nom "name" et composé des faces ayant les ids inclus dans "faces...
bool _readMetadata(IMesh *mesh, VtkFile &vtk_file)
Lecture des metadata.
void _readCellVariable(IMesh *mesh, VtkFile &vtk_file, const String &name_str, Integer nb_cell)
Permet de créer une variable aux mailles à partir des infos du fichier vtk.
bool _readStructuredGrid(IPrimaryMesh *mesh, VtkFile &, bool use_internal_partition)
Permet de lire un fichier vtk contenant une STRUCTURED_GRID.
void _readNodesUnstructuredGrid(IMesh *mesh, VtkFile &vtk_file, Array< Real3 > &node_coords)
Lecture des noeuds et de leur coordonnées.
void _readCellsUnstructuredGrid(IMesh *mesh, VtkFile &vtk_file, Array< Integer > &cells_nb_node, Array< ItemTypeId > &cells_type, Array< Int64 > &cells_connectivity)
Lecture des mailles et de leur connectivité.
void _readItemGroup(IMesh *mesh, VtkFile &vtk_file, const String &name_str, Integer nb_item, eItemKind ik, ConstArrayView< Int32 > local_id)
Permet de créer un groupe d'item.
void _readFacesMesh(IMesh *mesh, const String &file_name, const String &dir_name, bool use_internal_partition)
Permet de lire le fichier truc.vtkfaces.vtk (s'il existe).
bool _readData(IMesh *mesh, VtkFile &vtk_file, bool use_internal_partition, eItemKind cell_kind, Int32ConstArrayView local_id, Integer nb_node)
Permet de lire les données complémentaires (POINT_DATA / CELL_DATA).
bool readMesh(IPrimaryMesh *mesh, const String &file_name, const String &dir_name, bool use_internal_partition)
Permet de débuter la lecture d'un fichier vtk.
void _readNodeGroup(IMesh *mesh, VtkFile &vtk_file, const String &name, Integer nb_item)
Permet de créer un groupe de node.
bool _readUnstructuredGrid(IPrimaryMesh *mesh, VtkFile &vtk_file, bool use_internal_partition)
Permet de lire un fichier vtk contenant une UNSTRUCTURED_GRID.
eReturnType readMeshFromFile(IPrimaryMesh *mesh, const XmlNode &mesh_node, const String &file_name, const String &dir_name, bool use_internal_partition) override
Lit un maillage à partir d'un fichier.
bool allowExtension(const String &str) override
Vérifie si le service supporte les fichiers avec l'extension str.
Liste de noeuds d'un arbre DOM.
Integer size() const
Nombre d'éléments du vecteur.
Vue modifiable d'un tableau d'un type T.
void copy(Span< const T > rhs)
Copie les valeurs de rhs dans l'instance.
void resize(Int64 s)
Change le nombre d'éléments du tableau à s.
Vue constante d'un tableau de type T.
Interface du gestionnaire de traces.
virtual TraceMessage info()=0
Flot pour un message d'information.
Chaîne de caractères unicode.
ByteConstArrayView utf8() const
Retourne la conversion de l'instance dans l'encodage UTF-8.
const char * localstr() const
Retourne la conversion de l'instance dans l'encodage UTF-8.
Classe d'accès aux traces.
ITraceMng * traceMng() const
Gestionnaire de trace.
TraceMessage error() const
Flot pour un message d'erreur.
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flot pour un message de debug.
TraceMessage info() const
Flot pour un message d'information.
TraceMessage fatal() const
Flot pour un message d'erreur fatale.
Vecteur 1D de données avec sémantique par valeur (style STL).
#define ARCANE_REGISTER_SERVICE(aclass, a_service_property,...)
Macro pour enregistrer un service.
MeshVariableScalarRefT< Cell, Real > VariableCellReal
Grandeur au centre des mailles de type réel.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
@ ST_SubDomain
Le service s'utilise au niveau du sous-domaine.
eItemKind
Genre d'entité de maillage.
@ IK_Node
Entité de maillage de genre noeud.
@ IK_Cell
Entité de maillage de genre maille.
@ IK_Face
Entité de maillage de genre face.