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"
48#include "arcane/core/internal/VtkCellTypes.h"
60using namespace Arcane::VtkUtils;
119 VTK_MT_StructuredGrid,
120 VTK_MT_UnstructuredGrid
172 static const int BUFSIZE = 10000;
176 explicit VtkFile(std::istream* stream)
270 throw IOException(
"VtkFile::isEmptyNextLine()",
"Unexpected EndOfFile");
299 for (
int i = 0; i < BUFSIZE &&
m_buf[i] !=
'\0'; ++i) {
303 if (
m_buf[i] ==
'#') {
312 for (
int i = 0; i < BUFSIZE &&
m_buf[i] !=
'\0'; ++i) {
313 if (
m_buf[i] ==
'\r') {
329 throw IOException(
"VtkFile::isEmptyNextLine()",
"Not Good");
354 throw IOException(
"VtkFile::isEmptyNextLine()",
"Unexpected EndOfFile");
377 for (
int i = 0; i < BUFSIZE &&
m_buf[i] !=
'\0'; ++i) {
381 if (
m_buf[i] ==
'#') {
390 for (
int i = 0; i < BUFSIZE &&
m_buf[i] !=
'\0'; ++i) {
391 if (
m_buf[i] ==
'\r') {
399 throw IOException(
"VtkFile::getNextLine()",
"Not good");
418 (*m_stream) >> ws >> v;
423 throw IOException(
"VtkFile::getFloat()",
"Bad float");
442 (*m_stream) >> ws >> v;
447 throw IOException(
"VtkFile::getDouble()",
"Bad double");
466 (*m_stream) >> ws >> v;
486 constexpr size_t sizeofT =
sizeof(T);
496 for (
size_t i = 0; i <
sizeofT; i++) {
609 info() <<
"Titre du fichier VTK : " << title.
localstr();
645 error() <<
"Support exists only for 'STRUCTURED_GRID' and 'UNSTRUCTURED_GRID' formats (format=" <<
mesh_type_str <<
"')";
649 debug() <<
"Lecture en-tête OK";
653 case VTK_MT_StructuredGrid:
657 case VTK_MT_UnstructuredGrid:
659 debug() <<
"Lecture _readUnstructuredGrid OK";
663 debug() <<
"Lecture _readFacesMesh OK";
693 const char* buf =
nullptr;
699 std::istringstream
iline(buf);
704 error() <<
"Syntax error while reading grid dimensions";
721 std::istringstream
iline(buf);
726 error() <<
"Syntax error while reading grid dimensions";
766 for (Integer x = 0; x <
nb_node_x; ++x) {
767 for (Integer z = 0; z <
nb_node_z; ++z) {
768 for (Integer y = 0; y <
nb_node_y; ++y) {
800 for (Integer z = 0; z <
nb_cell_z; ++z) {
801 for (Integer y = 0; y <
nb_cell_y; ++y) {
802 for (Integer x = 0; x <
nb_cell_x; ++x) {
841 for (Integer z = 0; z <
nb_node_z; ++z) {
842 for (Integer y = 0; y <
nb_node_y; ++y) {
843 for (Integer x = 0; x <
nb_node_x; ++x) {
854 for (Integer z = 0; z <
nb_node_z; ++z) {
855 for (Integer y = 0; y <
nb_node_y; ++y) {
856 for (Integer x = 0; x <
nb_node_x; ++x) {
867 for (Integer z = 0; z <
nb_node_z; ++z) {
868 for (Integer y = 0; y <
nb_node_y; ++y) {
869 for (Integer x = 0; x <
nb_node_x; ++x) {
880 throw IOException(
"_readStructuredGrid",
"Invalid type name");
902 Integer face_local_id = face.
localId();
973 const char* func_name =
"VtkMeshIOService::_readNodesUnstructuredGrid()";
974 const char* buf =
vtk_file.getNextLine();
975 std::istringstream
iline(buf);
983 throw IOException(func_name,
"Syntax error while reading number of nodes");
988 throw IOException(A_FUNCINFO, String::format(
"Invalid number of nodes: n={0}",
nb_node));
996 for (Integer i = 0; i <
nb_node; ++i) {
1004 for (Integer i = 0; i <
nb_node; ++i) {
1012 for (Integer i = 0; i <
nb_node; ++i) {
1020 throw IOException(func_name,
"Invalid type name");
1046 ARCANE_UNUSED(mesh);
1048 const char* func_name =
"VtkMeshIOService::_readCellsUnstructuredGrid()";
1049 const char* buf =
vtk_file.getNextLine();
1054 std::istringstream
iline(buf);
1057 Integer nb_cell_node = 0;
1062 throw IOException(func_name,
"Syntax error while reading cells");
1066 if (
nb_cell < 0 || nb_cell_node < 0) {
1068 String::format(
"Invalid dimensions: nb_cell={0} nb_cell_node={1}",
1072 cells_nb_node.resize(
nb_cell);
1074 cells_connectivity.resize(nb_cell_node);
1078 for (Integer i = 0; i <
nb_cell; ++i) {
1080 cells_nb_node[i] = n;
1081 for (Integer
j = 0;
j < n; ++
j) {
1094 std::istringstream
iline(buf);
1100 throw IOException(func_name,
"Syntax error while reading cell types");
1110 for (Integer i = 0; i <
nb_cell; ++i) {
1112 Int16
it = vtkToArcaneCellType(
vtk_ct, cells_nb_node[i]);
1131 ARCANE_UNUSED(mesh);
1241 Integer nb_cell_node = 0;
1260 debug() <<
"Lecture _readNodesUnstructuredGrid OK";
1269 debug() <<
"Lecture _readCellsUnstructuredGrid OK";
1271 nb_cell = cells_nb_node.size();
1272 nb_cell_node = cells_connectivity.size();
1280 for (Integer i = 0; i <
nb_cell; ++i) {
1286 Int16
cell_dim =
itm->typeFromId(cells_type[i])->dimension();
1297 for (Int32 i = 0; i < 4; ++i)
1303 ARCANE_FATAL(
"The mesh contains cells of different dimension. nb0={0} nb1={1} nb2={2} nb3={3}",
1328 debug() <<
"Lecture _readData OK";
1350 std::ifstream
ifile(
file_name.localstr(), std::ifstream::binary);
1352 info() <<
"No face descriptor file found '" <<
file_name <<
"'";
1357 const char* buf = 0;
1361 info() <<
"Titre du fichier VTK : " << title;
1385 error() <<
"Face descriptor file type must be 'UNSTRUCTURED_GRID' (format=" <<
mesh_type_str <<
"')";
1468 pm->broadcast(
bb, 0);
1479 const char* buf = 0;
1485 debug() <<
"Read line";
1486 std::istringstream
iline(buf);
1523 for (Integer i = 0; i <
nb_fields; i++) {
1525 std::istringstream
iline(buf);
1537 if (
cstr ==
"GROUP_") {
1564 error() <<
"Expecting 'float' or 'double' data type, found=" <<
type_str;
1569 fatal() <<
"Unable to read POINT_DATA: feature not implemented";
1575 throw IOException(
"Unable to read face variables: feature not supported");
1598 if (
cstr ==
"GROUP_") {
1606 error() <<
"Expecting 'SCALARS' data type, found=" <<
data_str;
1631 error() <<
"Expecting 'float' or 'double' data type, found=" <<
type_str;
1638 fatal() <<
"Unable to read POINT_DATA: feature not implemented";
1644 throw IOException(
"Unable to read face variables: feature not supported");
1651 error() <<
"Expecting value CELL_DATA or POINT_DATA, found='" <<
data_str <<
"'";
1663 Integer len =
bv.size();
1671 String str = String::fromUtf8(bytes);
1672 info() <<
"FOUND STR=" << bytes.
size() <<
" " << str;
1681 info() <<
"Building variable: " << name;
1693 info() <<
"Building group: " << name;
1704 info() <<
"Create node group: " << name;
1726 info() <<
"Building face group '" << name <<
"'"
1749 mesh->
variableMng()->_internalApi()->addAutoDestroyVariable(
var);
1751 for (Integer i = 0; i <
nb_cell; ++i) {
1756 info() <<
"Variable build finished: " <<
vtk_file.isEof();
1777 info() <<
"Reading group info for group: " << name;
1780 for (Integer i = 0; i <
nb_item; ++i) {
1783 ids.add(local_id[i]);
1785 info() <<
"Building group: " << name <<
" nb_element=" << ids.size();
1807 info() <<
"Lecture infos groupes de noeuds pour le groupe: " << name;
1810 for (Integer i = 0; i <
nb_item; ++i) {
1815 info() <<
"Création groupe: " << name <<
" nb_element=" << ids.size();
1896 throw IOException(
"VtkMeshIOService::writeMeshToFile(): Unable to open file");
1897 ofile <<
"# vtk DataFile Version 2.0\n";
1898 ofile <<
"Maillage Arcane\n";
1900 ofile <<
"DATASET UNSTRUCTURED_GRID\n";
1933 Integer item_nb_node = item.
nbNode();
1934 ofile << item_nb_node;
1945 ofile << type <<
'\n';
1963void VtkLegacyMeshWriter::
1966 info() <<
"Saving groups for family name=" << family->
name();
1974 if (group.
name() ==
"OuterFaces")
1976 ofile <<
"SCALARS GROUP_" << group.
name() <<
" int 1\n";
1977 ofile <<
"LOOKUP_TABLE default\n";
2062 m_trace_mng->
info() <<
"VtkLegacy Reader (ICaseMeshReader) file_name=" <<
fname;
2063 bool ret =
vtk_service.readMesh(pm,
fname, m_read_info.directoryName(), m_read_info.isParallelRead());
2087 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.