NEMoSys  0.63.0
A modular, extensible resource with robust automated mesh generation, mesh quality analysis, adaptive mesh refinement, and data transfer between arbitrary meshes.
vtkMesh.H
Go to the documentation of this file.
1 /*******************************************************************************
2 * Promesh *
3 * Copyright (C) 2022, IllinoisRocstar LLC. All rights reserved. *
4 * *
5 * Promesh is the property of IllinoisRocstar LLC. *
6 * *
7 * IllinoisRocstar LLC *
8 * Champaign, IL *
9 * www.illinoisrocstar.com *
10 * promesh@illinoisrocstar.com *
11 *******************************************************************************/
12 /*******************************************************************************
13 * This file is part of Promesh *
14 * *
15 * This version of Promesh is free software: you can redistribute it and/or *
16 * modify it under the terms of the GNU Lesser General Public License as *
17 * published by the Free Software Foundation, either version 3 of the License, *
18 * or (at your option) any later version. *
19 * *
20 * Promesh is distributed in the hope that it will be useful, but WITHOUT ANY *
21 * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS *
22 * FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more *
23 * details. *
24 * *
25 * You should have received a copy of the GNU Lesser General Public License *
26 * along with this program. If not, see <https://www.gnu.org/licenses/>. *
27 * *
28 *******************************************************************************/
29 #ifndef NEMOSYS_VTKMESH_H_
30 #define NEMOSYS_VTKMESH_H_
31 
32 #include "nemosys_export.h"
33 #include "Mesh/meshBase.H"
34 
35 #include <vtkUnstructuredGrid.h>
36 
37 #ifndef SWIG
38 // read files from any vtk format (except legacy .vtk)
39 template <class TReader>
40 vtkDataSet *ReadAnXMLOrSTLFile(const std::string &fileName) {
41  vtkSmartPointer<TReader> reader = vtkSmartPointer<TReader>::New();
42  reader->SetFileName(fileName.c_str());
43  reader->Update();
44  reader->GetOutput()->Register(reader);
45  return vtkDataSet::SafeDownCast(reader->GetOutput());
46 }
47 
48 // read an ASCII legacy vtk file containing an unstructured grid
49 // the available vtk readers do not suffice for this purpose if the file
50 // contains more than one point or cell data array. This function has mainly
51 // been tested on vtk output from MFEM
52 vtkSmartPointer<vtkUnstructuredGrid> ReadALegacyVTKFile(
53  const std::string &fileName);
54 
55 vtkSmartPointer<vtkUnstructuredGrid> ReadDegenerateVTKFile(
56  const std::string &fileName);
57 
58 // helpers for reading legacy vtk
59 
60 // helper that reads/checks header from legacy vtk
61 bool readLegacyVTKHeader(const std::string &line);
62 // helper that reads field data from legacy vtk
63 bool readLegacyVTKFieldData(const std::istream &meshStream,
64  const std::string &line,
65  vtkSmartPointer<vtkUnstructuredGrid> dataSet_tmp);
66 // helper that reads points from legacy vtk
67 bool readLegacyVTKPoints(const std::istream &meshStream,
68  const std::string &line, nemId_t &numPoints,
69  vtkSmartPointer<vtkPoints> points,
70  vtkSmartPointer<vtkUnstructuredGrid> dataSet_tmp);
71 // helper that reads cells from legacy vtk
73  const std::istream &meshStream, const std::string &line, nemId_t &numCells,
74  const std::vector<vtkSmartPointer<vtkIdList>> &vtkCellIds,
75  vtkSmartPointer<vtkUnstructuredGrid> dataSet_tmp);
76 // helper that reads point and cell data from legacy vtk
77 bool readLegacyVTKData(const std::ifstream &meshStream, const std::string &line,
78  nemId_t numTuple, bool pointOrCell, bool &hasPointOrCell,
79  vtkSmartPointer<vtkUnstructuredGrid> dataSet_tmp);
80 // helper that casts data arrays to type specified in legacy file (not used rn)
81 void addLegacyVTKData(vtkDataArray *arr, const std::string &type,
82  bool pointOrCell,
83  vtkSmartPointer<vtkUnstructuredGrid> dataSet_tmp);
84 
85 // method to write vtk grids
86 template <class TWriter>
87 void writeVTFile(const std::string &fname,
88  vtkSmartPointer<vtkDataSet> dataSet) {
89  vtkSmartPointer<TWriter> Writer = vtkSmartPointer<TWriter>::New();
90  Writer->SetFileName(fname.c_str());
91  Writer->SetInputData(dataSet);
92  // WARNING: ASCII does not work for STL and legacy VTK files
93  // WARNING: large files will be generated
94  // Writer->SetDataModeToAscii();
95  Writer->Update();
96  Writer->Write();
97 }
98 #endif // SWIG
99 
100 class NEMOSYS_DEPRECATED_EXPORT vtkMesh
101  : public meshBase {
102  // constructor and destructor
103  public:
104  vtkMesh() = default;
105  explicit vtkMesh(const std::string &fname);
106  // limited support for conversion between formats
107  vtkMesh(const std::string &fname1, const std::string &fname2);
108  // put existing vtkDataSet into meshBase
109  vtkMesh(vtkSmartPointer<vtkDataSet> dataSet_tmp, const std::string &fname);
110  /* create from coordinates and connectivities.
111  NOTE: use of this is only valid when mesh has one cell type.
112  cellType one of the vtkCellType enums.
113  Currently, only VTK_TETRA and VTK_TRIANGLE are supported */
114  vtkMesh(const std::vector<double> &xCrds, const std::vector<double> &yCrds,
115  const std::vector<double> &zCrds,
116  const std::vector<nemId_t> &elemConn, int cellType,
117  const std::string &newname);
118 
119  ~vtkMesh() override { std::cout << "vtkMesh destroyed" << std::endl; }
120 
121  // access
122  public:
123  // abstract read method reserved for derived classes
124  void read(const std::string &fname) override {}
125  // merge vtk mesh with current mesh
126  void merge(vtkSmartPointer<vtkDataSet> dataSet_new);
127  // get point with id
128  std::vector<double> getPoint(nemId_t id) const override;
129  // get 3 vecs with x,y and z coords
130  std::vector<std::vector<double>> getVertCrds() const override;
131  // get cell with id : returns point indices and respective coordinates
132  std::map<nemId_t, std::vector<double>> getCell(nemId_t id) const override;
133  std::vector<std::vector<double>> getCellVec(nemId_t id) const override;
134  // get diameter of circumsphere of each cell
135  std::vector<double> getCellLengths() const override;
136  // get center of a cell
137  std::vector<double> getCellCenter(nemId_t cellID) const override;
138  // get cell type as an integer
139  // assumes all elements are the same type
140  int getCellType() const override;
141  // get edge lengths of dataSet
142  void inspectEdges(const std::string &ofname) const override;
143  std::vector<nemId_t> getConnectivities() const override;
144 
145  // processing
146  public:
147  vtkSmartPointer<vtkDataSet> extractSurface() override;
148 
149  // diagnostics
150  public:
151  void report() const override;
152  void write() const override { meshBase::write(); }
153  void write(const std::string &fname) const override;
154 
155  // set and get point and cell data
156  public:
157  // set point data (numComponents per point determined by dim of data[0]
158  void setPointDataArray(const std::string &name,
159  const std::vector<std::vector<double>> &data) override;
160  void setPointDataArray(const std::string &name,
161  const std::vector<double> &data) override;
162  // set cell data (numComponents per cell determined by dim of data[0])
163  void setCellDataArray(const std::string &name,
164  const std::vector<std::vector<double>> &data) override;
165  // set scalar cell data
166  void setCellDataArray(const std::string &name,
167  const std::vector<double> &data) override;
168  // remove point data with given id from dataSet if it exists
169  void unsetPointDataArray(int arrayID) override;
170  void unsetPointDataArray(const std::string &name) override;
171  // remove cell data with given id from dataSet if it exists
172  void unsetCellDataArray(int arrayID) override;
173  void unsetCellDataArray(const std::string &name) override;
174  // remove field data with given id from dataSet
175  void unsetFieldDataArray(const std::string &name) override;
176  // get scalar point or cell data array. assumes data is not allocated prior to
177  // calling
178  void getPointDataArray(const std::string &name, std::vector<double> &data)
179  override;
180  void getPointDataArray(int arrayId, std::vector<double> &data) override;
181  int getCellDataIdx(const std::string &name) override;
182  void getCellDataArray(const std::string &name, std::vector<double> &data)
183  override;
184  void getCellDataArray(int arrayId, std::vector<double> &data) override;
185 };
186 
187 #endif // NEMOSYS_VTKMESH_H_
vtkSmartPointer< vtkUnstructuredGrid > ReadALegacyVTKFile(const std::string &fileName)
Definition: vtkMesh.C:502
bool readLegacyVTKHeader(const std::string &line)
Definition: vtkMesh.C:294
data_type data
Edge/face with sorted point ids (a, b, c, ...) is located at some index i in data[b], with data[b][i].first == [a, c] (for edges, third point id treated as -1).
virtual void inspectEdges(const std::string &ofname) const =0
get edge lengths of dataSet
virtual void getPointDataArray(const std::string &name, std::vector< double > &data)
get scalar point or cell data array.
Definition: meshBase.H:349
virtual void setPointDataArray(const std::string &name, const std::vector< std::vector< double >> &data)
register data to dataSet&#39;s point data
Definition: meshBase.H:319
void addLegacyVTKData(vtkDataArray *arr, const std::string &type, bool pointOrCell, vtkSmartPointer< vtkUnstructuredGrid > dataSet_tmp)
virtual std::vector< double > getCellLengths() const =0
get diameter of circumsphere of each cell
vtkSmartPointer< vtkUnstructuredGrid > ReadDegenerateVTKFile(const std::string &fileName)
Definition: vtkMesh.C:581
A brief description of meshBase.
Definition: meshBase.H:64
std::size_t nemId_t
Definition: meshBase.H:51
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
virtual void setCellDataArray(const std::string &name, const std::vector< std::vector< double >> &data)
register data to dataSet&#39;s cell data
Definition: meshBase.H:334
void writeVTFile(const std::string &fname, vtkSmartPointer< vtkDataSet > dataSet)
Definition: vtkMesh.H:87
~vtkMesh() override
Definition: vtkMesh.H:119
virtual int getCellDataIdx(const std::string &name)
<>
Definition: meshBase.H:363
virtual std::vector< std::vector< double > > getCellVec(nemId_t id) const =0
get vector of coords of cell with id
bool readLegacyVTKCells(const std::istream &meshStream, const std::string &line, nemId_t &numCells, const std::vector< vtkSmartPointer< vtkIdList >> &vtkCellIds, vtkSmartPointer< vtkUnstructuredGrid > dataSet_tmp)
virtual int getCellType() const =0
get cell type as an integer assumes all elements are the same type
virtual std::vector< nemId_t > getConnectivities() const =0
get connectivities.
virtual void unsetPointDataArray(int arrayID)
delete array with id from dataSet&#39;s point data
Definition: meshBase.H:381
virtual std::map< nemId_t, std::vector< double > > getCell(nemId_t id) const =0
get cell with id
virtual vtkSmartPointer< vtkDataSet > extractSurface()=0
extract the surface mesh
virtual void write() const
write the mesh to file named after the private var &#39;filename&#39;.
Definition: meshBase.H:598
VTKCellType cellType
Definition: inpGeoMesh.C:129
bool readLegacyVTKPoints(const std::istream &meshStream, const std::string &line, nemId_t &numPoints, vtkSmartPointer< vtkPoints > points, vtkSmartPointer< vtkUnstructuredGrid > dataSet_tmp)
void read(const std::string &fname) override
abstract read method reserved for derived classes
Definition: vtkMesh.H:124
bool readLegacyVTKFieldData(const std::istream &meshStream, const std::string &line, vtkSmartPointer< vtkUnstructuredGrid > dataSet_tmp)
virtual void unsetCellDataArray(int arrayID)
delete array with id from dataSet&#39;s cell data
Definition: meshBase.H:391
std::vector< vtkIdType > points
points given by id in .inp file
Definition: inpGeoMesh.C:133
vtkDataSet * ReadAnXMLOrSTLFile(const std::string &fileName)
Definition: vtkMesh.H:40
virtual std::vector< std::vector< double > > getVertCrds() const =0
get 3 vecs with x,y and z coords
virtual std::vector< double > getCellCenter(nemId_t cellID) const =0
get center of a cell
virtual void report() const
generate a report of the mesh
Definition: meshBase.H:540
virtual void getCellDataArray(const std::string &name, std::vector< double > &data)
<>
Definition: meshBase.H:369
bool readLegacyVTKData(const std::ifstream &meshStream, const std::string &line, nemId_t numTuple, bool pointOrCell, bool &hasPointOrCell, vtkSmartPointer< vtkUnstructuredGrid > dataSet_tmp)
void write() const override
write the mesh to file named after the private var &#39;filename&#39;.
Definition: vtkMesh.H:152
virtual std::vector< double > getPoint(nemId_t id) const =0
get point with id
virtual void unsetFieldDataArray(const std::string &name)
delete array with id from dataSet&#39;s field data
Definition: meshBase.H:401