1 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 // + This file is part of enGrid. +
5 // + Copyright 2008-2014 enGits GmbH +
7 // + enGrid is free software: you can redistribute it and/or modify +
8 // + it under the terms of the GNU General Public License as published by +
9 // + the Free Software Foundation, either version 3 of the License, or +
10 // + (at your option) any later version. +
12 // + enGrid is distributed in the hope that it will be useful, +
13 // + but WITHOUT ANY WARRANTY; without even the implied warranty of +
14 // + MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the +
15 // + GNU General Public License for more details. +
17 // + You should have received a copy of the GNU General Public License +
18 // + along with enGrid. If not, see <http://www.gnu.org/licenses/>. +
20 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
21 #include "vtkreader.h"
23 #include <vtkUnstructuredGridReader.h>
26 #include "guimainwindow.h"
28 VtkReader::VtkReader()
30 setFormat("legacy VTK files(*.vtk)");
34 void VtkReader::createBoundaryFaces()
36 QList
<QVector
<vtkIdType
> > all_faces
;
37 MeshPartition
part(m_Grid
, true);
38 for (vtkIdType id_cell
= 0; id_cell
< m_Grid
->GetNumberOfCells(); ++id_cell
) {
39 for (int i
= 0; i
< part
.c2cGSize(id_cell
); ++i
) {
40 if (part
.c2cGG(id_cell
, i
) < 0) {
41 QVector
<vtkIdType
> face
;
42 getFaceOfCell(m_Grid
, id_cell
, i
, face
);
47 EG_VTKSP(vtkUnstructuredGrid
, grid
);
48 allocateGrid(grid
, m_Grid
->GetNumberOfCells() + all_faces
.size(), m_Grid
->GetNumberOfPoints());
49 for (vtkIdType id_node
= 0; id_node
< m_Grid
->GetNumberOfPoints(); ++id_node
) {
51 m_Grid
->GetPoint(id_node
, x
.data());
52 grid
->GetPoints()->SetPoint(id_node
, x
.data());
54 for (vtkIdType id_cell
= 0; id_cell
< m_Grid
->GetNumberOfCells(); ++id_cell
) {
55 copyCell(m_Grid
, id_cell
, grid
);
57 EG_VTKDCC(vtkIntArray
, bc
, grid
, "cell_code");
58 foreach (QVector
<vtkIdType
> face
, all_faces
) {
59 EG_VTKSP(vtkIdList
, pts
);
60 pts
->SetNumberOfIds(face
.size());
61 for (int i
= 0; i
< face
.size(); ++i
) {
62 pts
->SetId(i
, face
[i
]);
65 if (face
.size() == 3) {
66 id_cell
= grid
->InsertNextCell(VTK_TRIANGLE
, pts
);
67 } else if (face
.size() == 4) {
68 id_cell
= grid
->InsertNextCell(VTK_QUAD
, pts
);
70 id_cell
= grid
->InsertNextCell(VTK_POLYGON
, pts
);
72 bc
->SetValue(id_cell
, 1);
74 makeCopy(grid
, m_Grid
);
75 GuiMainWindow::pointer()->addBC(1, BoundaryCondition("all_faces", "patch"));
78 void VtkReader::operate()
81 QFileInfo
file_info(GuiMainWindow::pointer()->getFilename());
82 readInputFileName(file_info
.completeBaseName() + ".vtk");
84 EG_VTKSP(vtkUnstructuredGridReader
,vtk
);
85 vtk
->SetFileName(qPrintable(getFileName()));
87 //m_Grid->DeepCopy(vtk->GetOutput());
88 makeCopy(vtk
->GetOutput(), m_Grid
);
89 createBoundaryFaces();
90 createBasicFields(m_Grid
, m_Grid
->GetNumberOfCells(), m_Grid
->GetNumberOfPoints());
91 UpdateNodeIndex(m_Grid
);
92 UpdateCellIndex(m_Grid
);