CEPS  24.01
Cardiac ElectroPhysiology Simulator
VtkMeshReader.hpp
Go to the documentation of this file.
1 /*~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2  This file is part of CEPS.
3 
4  CEPS is free software: you can redistribute it and/or modify
5  it under the terms of the GNU General Public License as published by
6  the Free Software Foundation, either version 3 of the License, or
7  (at your option) any later version.
8 
9  CEPS is distributed in the hope that it will be useful,
10  but WITHOUT ANY WARRANTY; without even the implied warranty of
11  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  GNU General Public License for more details.
13 
14  You should have received a copy of the GNU General Public License
15  along with CEPS (see file LICENSE at root of project).
16  If not, see <https://www.gnu.org/licenses/>.
17 
18 
19  Copyright 2019-2024 Inria, Universite de Bordeaux
20 
21  Authors, in alphabetical order:
22 
23  Pierre-Elliott BECUE, Florian CARO, Yves COUDIERE(*), Andjela DAVIDOVIC,
24  Charlie DOUANLA-LONTSI, Marc FUENTES, Mehdi JUHOOR, Michael LEGUEBE(*),
25  Pauline MIGERDITICHAN, Valentin PANNETIER(*), Nejib ZEMZEMI.
26  * : currently active authors
27 
28 ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~*/
30 #pragma once
31 
33 
34 #include <vtkCellData.h>
35 #include <vtkDataSet.h>
36 #include <vtkDataSetReader.h>
37 #include <vtkGenericCell.h>
38 #include <vtkIdList.h>
39 #include <vtkPointData.h>
40 #include <vtkPointSet.h>
41 #include <vtkPoints.h>
42 #include <vtkPolyData.h>
43 #include <vtkSmartPointer.h>
44 #include <vtkUnstructuredGrid.h>
45 #include <vtkXMLPUnstructuredGridReader.h>
46 #include <vtkCommand.h>
47 #include <vtkNew.h>
48 
57 class VtkMeshReader : public MeshReader
58 {
59 
60  public:
61 
63  VtkMeshReader(const CepsString &fileName, const CepsUInt &dim, CepsBool quads = false);
64 
66  ~VtkMeshReader() override;
67 
70  void
71  initialize() override;
72 
74  CepsBool
75  isVtkFile(const CepsString &fileName) const;
76 
78  void
79  readNodes(
80  CepsUInt indexStart,
81  CepsUInt indexEnd,
82  CepsUInt nodeOffset,
83  CepsVector<CepsReal>& coords,
85  CepsVector<CepsUInt>& attrPtr
86  ) override;
87 
89  void
90  readCells(
91  CepsUInt indexStart,
92  CepsUInt indexEnd,
93  CepsUInt nodeOffset,
94  CepsUInt cellOffset,
96  CepsVector<CepsUInt>& nodeIndicesPtr,
97  CepsVector<CepsAttribute>& cellAttr,
98  CepsVector<CepsUInt>& cellAttrPtr,
99  CepsVector<CepsChar>& isBoundaryCell
100  ) override;
101 
102 
103  private:
104 
106  void
107  getNextCell(
108  CepsUInt cellId,
109  CepsUInt fileCellId,
110  CepsUInt nodeOffset,
111  CepsVector<CepsNodeGlobalIndex>& nodeIndices,
112  CepsVector<CepsUInt>& nodeIndicesPtr,
113  CepsVector<CepsAttribute>& cellAttr,
114  CepsVector<CepsUInt>& cellAttrPtr,
115  CepsVector<CepsChar>& isBoundaryCell,
116  CepsBool boundaryFlag
117  ) override;
118 
120  void
121  initializePvtu(const CepsString &fileName);
122 
124  void
125  initializeLegacy(const CepsString &fileName);
126 
128  void
130 
132  void
134 
135 
136  protected:
137 
138  vtkSmartPointer<vtkDataSet> m_mesh;
139 
142 };
std::basic_string< CepsChar > CepsString
C++ format string.
Definition: CepsTypes.hpp:128
std::vector< _Type, _Alloc > CepsVector
C++ vector.
Definition: CepsTypes.hpp:155
bool CepsBool
Booleans.
Definition: CepsTypes.hpp:124
std::make_unsigned_t< CepsInt > CepsUInt
Unsigned version on CepsInt.
Definition: CepsTypes.hpp:109
Abstract base class that encapsulates primary functionalities of each mesh reader.
Definition: MeshReader.hpp:42
This class enables reading of vtk dataset files.
~VtkMeshReader() override
Destructor.
void getNextCell(CepsUInt cellId, CepsUInt fileCellId, CepsUInt nodeOffset, CepsVector< CepsNodeGlobalIndex > &nodeIndices, CepsVector< CepsUInt > &nodeIndicesPtr, CepsVector< CepsAttribute > &cellAttr, CepsVector< CepsUInt > &cellAttrPtr, CepsVector< CepsChar > &isBoundaryCell, CepsBool boundaryFlag) override
Reads a single cell in mesh file.
void readNodes(CepsUInt indexStart, CepsUInt indexEnd, CepsUInt nodeOffset, CepsVector< CepsReal > &coords, CepsVector< CepsAttribute > &attr, CepsVector< CepsUInt > &attrPtr) override
Reads several nodes in mesh file.
void determinePointDataToRead()
Selects the integer array that should be present, which contains node attributes.
CepsVector< vtkDataArray * > m_cellArrayToRead
Arrays of data defined on cells.
void readCells(CepsUInt indexStart, CepsUInt indexEnd, CepsUInt nodeOffset, CepsUInt cellOffset, CepsVector< CepsNodeGlobalIndex > &nodeIndices, CepsVector< CepsUInt > &nodeIndicesPtr, CepsVector< CepsAttribute > &cellAttr, CepsVector< CepsUInt > &cellAttrPtr, CepsVector< CepsChar > &isBoundaryCell) override
Reads several cells in mesh file.
VtkMeshReader(const CepsString &fileName, const CepsUInt &dim, CepsBool quads=false)
Constructor with full file name.
void determineCellDataToRead()
Selects the integer array that should be present, which contains cell attributes.
void initializePvtu(const CepsString &fileName)
Reads global information for parallel unstructured grids.
vtkSmartPointer< vtkDataSet > m_mesh
VTK mesh structure.
CepsVector< vtkDataArray * > m_pointArrayToRead
Arrays of data defined on points.
void initializeLegacy(const CepsString &fileName)
Reads global information for legacy files.
void initialize() override
Read nodes from indexStart to indexEnd, filling the arrays of coordinates and attributes.
CepsBool isVtkFile(const CepsString &fileName) const
Test if file is legacy file .vtk or .pvtu.