Skip to content
Snippets Groups Projects
Select Git revision
  • 06b44fdf03912eedd32c0642b8f92f320d7a1953
  • develop default protected
  • feature/advection
  • feature/composite-scheme-other-fluxes
  • origin/stage/bouguettaia
  • save_clemence
  • feature/local-dt-fsi
  • feature/variational-hydro
  • feature/gmsh-reader
  • feature/reconstruction
  • feature/kinetic-schemes
  • feature/composite-scheme-sources
  • feature/serraille
  • feature/composite-scheme
  • hyperplastic
  • feature/polynomials
  • feature/gks
  • feature/implicit-solver-o2
  • feature/coupling_module
  • feature/implicit-solver
  • feature/merge-local-dt-fsi
  • v0.5.0 protected
  • v0.4.1 protected
  • v0.4.0 protected
  • v0.3.0 protected
  • v0.2.0 protected
  • v0.1.0 protected
  • Kidder
  • v0.0.4 protected
  • v0.0.3 protected
  • v0.0.2 protected
  • v0 protected
  • v0.0.1 protected
33 results

GmshReader.hpp

Blame
  • Stephane Del Pino's avatar
    Stéphane Del Pino authored
    'typedef ConcreteType AliasType;' -> 'using AliasType = ConcreteType;'
    
    This is much cleaner
    06b44fdf
    History
    GmshReader.hpp 4.82 KiB
    #ifndef GMSH_READER_HPP
    #define GMSH_READER_HPP
    
    #include <Array.hpp>
    #include <TinyVector.hpp>
    
    #include <Mesh.hpp>
    
    #include <RefId.hpp>
    
    #include <string>
    #include <array>
    #include <vector>
    #include <map>
    #include <fstream>
    
    class GmshReader
    {
    public:
      using R3 = TinyVector<3, double>;
    
      class PhysicalRefId
      {
       private:
        int m_dimension;
        RefId m_ref_id;
    
       public:
        friend std::ostream& operator<<(std::ostream& os,
                                        const PhysicalRefId& physical_ref_id)
        {
          os << physical_ref_id.m_ref_id << "[dimension " << physical_ref_id.m_dimension << "]";
          return os;
        }
    
        bool operator<(const PhysicalRefId& physical_ref_id) const
        {
          return m_ref_id.tagNumber()<physical_ref_id.m_ref_id.tagNumber();
        }
    
        bool operator==(const PhysicalRefId& physical_ref_id) const
        {
          return m_ref_id.tagNumber()==physical_ref_id.m_ref_id.tagNumber();
        }
    
        const int& dimension() const
        {
          return m_dimension;
        }
    
        const RefId& refId() const
        {
          return m_ref_id;
        }
    
        PhysicalRefId& operator=(const PhysicalRefId&) = default;
        PhysicalRefId& operator=(PhysicalRefId&&) =default;
        PhysicalRefId(const int& dimension,
                      const RefId& ref_id)
            : m_dimension(dimension),
              m_ref_id(ref_id)
        {
          ;
        }
        PhysicalRefId() = default;
        PhysicalRefId(const PhysicalRefId&) = default;
        PhysicalRefId(PhysicalRefId&&) = default;
        ~PhysicalRefId() = default;
      };
    
    private:
      std::ifstream m_fin;
      const std::string m_filename;
    
      /**
       * Gmsh format provides a numbered, none ordrered and none dense
       * vertices list, this stores the number of read vertices.
       */
      std::vector<int> __verticesNumbers;
    
      Array<R3> __vertices;
    
      using Point = unsigned int;
      Array<Point> __points;
      std::vector<int> __points_ref;
    
      using Edge = TinyVector<2,unsigned int>;
      Array<Edge> __edges;
      std::vector<int> __edges_ref;
    
      using Triangle = TinyVector<3,unsigned int>;
      Array<Triangle> __triangles;
      std::vector<int> __triangles_ref;
    
      using Quadrangle = TinyVector<4,unsigned int>;
      Array<Quadrangle> __quadrangles;
      std::vector<int> __quadrangles_ref;
    
      using Tetrahedron = TinyVector<4,unsigned int>;
      Array<Tetrahedron> __tetrahedra;
      std::vector<int> __tetrahedra_ref;
    
      using Hexahedron = TinyVector<8,unsigned int>;
      Array<Hexahedron> __hexahedra;
      std::vector<int> __hexahedra_ref;
    
      /**
       * Gmsh format provides a numbered, none ordrered and none dense
       * vertices list, this provides vertices renumbering correspondance
       */
      std::vector<int> __verticesCorrepondance;
    
      /**
       * elements types
       */
      std::vector<short> __elementType;
    
      /**
       * References
       */
      std::vector<int> __references;
    
      /**
       * References
       */
      std::vector<std::vector<int> > __elementVertices;
    
      /**
       * Stores the number of nodes associated to each primitive
       *
       */
      std::vector<int> __numberOfPrimitiveNodes;
    
      /**
       * Array of boolean describing gmsh supported primitives by ff3d
       *
       */
      std::array<bool, 15> __supportedPrimitives;
    
      /**
       * Primitives names according to there number
       *
       */
      std::map<int, std::string> __primitivesNames;
    
      bool __binary;		/**< true if binary format */
      bool __convertEndian;		/**< true if needs to adapt endianess */
    
      std::map<unsigned int, PhysicalRefId> m_physical_ref_map;
    
      int _getInteger()
      {
        int i;
        m_fin >> i;
    
        return std::move(i);
      }
    
      double _getReal()
      {
        double d;
        m_fin >> d;
    
        return std::move(d);
      }
    
      /**
       * Adapts gmsh data to ff3d's structures
       *
       */
      void __proceedData();
    
      /**
       * Reads data in format 2.2
       *
       */
      void __readGmshFormat2_2();
    
      /**
       * List of allowed keyword in mesh file
       *
       */
      enum KeywordType {
        // gmsh 2.2
        MESHFORMAT,
        ENDMESHFORMAT,
        NODES,
        ENDNODES,
        ELEMENTS,
        ENDELEMENTS,
        PHYSICALNAMES,
        ENDPHYSICALNAMES,
    
        Unknown,
        EndOfFile
      };
    
      /**
       * List of known keywords
       *
       */
      using KeywordList = std::map<std::string, int>;
    
      KeywordList __keywordList;	/**< The keyword list */
    
      /**
       * Type for keyword
       *
       */
      using Keyword = std::pair<std::string, int>;
    
      /**
       * Skips next comments if exists
       *
       */
      void __skipComments();
    
      /**
       * Reads the next keyword and returns its token
       *
       * @return KeywordToken
       */
      Keyword __nextKeyword();
    
      /**
       * get list of vertices
       *
       */
      void __getVertices();
    
      /**
       * Read list of vertices
       *
       */
      void __readVertices();
    
      /**
       * Read all elements in format 2.0
       *
       */
      void __readElements2_2();
    
      /**
       * Reads physical names
       *
       */
      void __readPhysicalNames2_2();
    
      std::shared_ptr<IMesh> m_mesh;
    public:
    
      std::shared_ptr<IMesh> mesh() const
      {
        return m_mesh;
      }
    
      GmshReader(const std::string& filename);
      ~GmshReader() = default;
    
    
    };
    
    #endif // GMSH_READER_HPP