deactivated deprecated start_engrid copy
[engrid.git] / src / blenderreader.cpp
blob5655f9a1602230964df995c0d2d29e3586a6a42c
1 //
2 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 // + +
4 // + This file is part of enGrid. +
5 // + +
6 // + Copyright 2008,2009 Oliver Gloth +
7 // + +
8 // + enGrid is free software: you can redistribute it and/or modify +
9 // + it under the terms of the GNU General Public License as published by +
10 // + the Free Software Foundation, either version 3 of the License, or +
11 // + (at your option) any later version. +
12 // + +
13 // + enGrid is distributed in the hope that it will be useful, +
14 // + but WITHOUT ANY WARRANTY; without even the implied warranty of +
15 // + MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the +
16 // + GNU General Public License for more details. +
17 // + +
18 // + You should have received a copy of the GNU General Public License +
19 // + along with enGrid. If not, see <http://www.gnu.org/licenses/>. +
20 // + +
21 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
24 #include "blenderreader.h"
25 #include "guimainwindow.h"
27 BlenderReader::BlenderReader()
29 setFormat("Blender/Engrid files(*.begc *.BEGC)");
32 void BlenderReader::operate()
34 try {
35 QFileInfo file_info(GuiMainWindow::pointer()->getFilename());
36 readInputFileName(file_info.completeBaseName() + ".begc");
37 if (isValid()) {
38 // read raw data from exported file
39 QFile file(getFileName());
40 file.open(QIODevice::ReadOnly | QIODevice::Text);
41 QTextStream f(&file);
42 QList<vec3_t> rnodes;
43 QList<QVector<int> > rfaces;
44 int num_parts;
45 f >> num_parts;
46 QVector<QString> part_name(num_parts);
47 for (int i_part = 0; i_part < num_parts; ++i_part) {
48 f >> part_name[i_part];
50 QVector<QString> sorted_part_name = part_name;
51 qSort(sorted_part_name);
52 QVector<int> part_bc(part_name.size());
53 for (int i_part = 0; i_part < num_parts; ++i_part) {
54 part_bc[i_part] = sorted_part_name.indexOf(part_name[i_part]) + 1;
56 for (int i_part = 0; i_part < num_parts; ++i_part) {
57 int num_nodes, num_faces;
58 f >> num_nodes >> num_faces;
59 for (int i = 0; i < num_nodes; ++i) {
60 vec3_t x;
61 f >> x[0] >> x[1] >> x[2];
62 rnodes.push_back(x);
64 for (int i = 0; i < num_faces; ++i) {
65 int N;
66 f >> N;
67 QVector<int> face(N+1);
68 face[0] = i_part;
69 for (int j = 0; j < N; ++j) {
70 f >> face[j+1];
72 rfaces.push_back(face);
75 QVector<vec3_t> nodes(rnodes.size());
76 qCopy(rnodes.begin(), rnodes.end(), nodes.begin());
77 QVector<QVector<int> > faces(rfaces.size());
78 qCopy(rfaces.begin(), rfaces.end(), faces.begin());
80 // find smallest edge length
81 double L = 1e99;
82 foreach (QVector<int> face, faces) {
83 for (int i = 1; i < face.size(); ++i) {
84 int n1 = face[i];
85 int n2 = face[1];
86 if (i < face.size() - 1) {
87 n2 = face[i+1];
89 double l = (nodes[n1] - nodes[n2]).abs();
90 L = min(l, L);
94 // delete duplicate nodes
95 QList<vec3_t> non_dup;
96 QVector<int> o2n(nodes.size());
97 int num_non_dup = 0;
98 for (int i = 0; i < nodes.size(); ++i) {
99 o2n[i] = num_non_dup;
100 bool dup = false;
101 for (int j = 0; j < i; ++j) {
102 if ((nodes[i] - nodes[j]).abs() < 0.1*L) {
103 o2n[i] = o2n[j];
104 dup = true;
105 break;
108 if (!dup) {
109 non_dup.push_back(nodes[i]);
110 ++num_non_dup;
114 allocateGrid(m_Grid, faces.size(), non_dup.size());
115 EG_VTKDCC(vtkIntArray, cell_code, m_Grid, "cell_code");
116 EG_VTKDCC(vtkIntArray, orgdir, m_Grid, "cell_orgdir");
117 EG_VTKDCC(vtkIntArray, voldir, m_Grid, "cell_voldir");
118 EG_VTKDCC(vtkIntArray, curdir, m_Grid, "cell_curdir");
119 vtkIdType id_node = 0;
120 foreach (vec3_t x, non_dup) {
121 m_Grid->GetPoints()->SetPoint(id_node, x.data());
122 ++id_node;
125 foreach (QVector<int> face, faces) {
126 if (face.size() == 4) {
127 vtkIdType pts[3];
128 pts[0] = o2n[face[1]];
129 pts[1] = o2n[face[2]];
130 pts[2] = o2n[face[3]];
131 vtkIdType id_cell = m_Grid->InsertNextCell(VTK_TRIANGLE, 3, pts);
132 cell_code->SetValue(id_cell, part_bc[face[0]]);
133 orgdir->SetValue(id_cell, 0);
134 voldir->SetValue(id_cell, 0);
135 curdir->SetValue(id_cell, 0);
137 if (face.size() == 5) {
138 vtkIdType pts[4];
139 pts[0] = o2n[face[1]];
140 pts[1] = o2n[face[2]];
141 pts[2] = o2n[face[3]];
142 pts[3] = o2n[face[4]];
143 vtkIdType id_cell = m_Grid->InsertNextCell(VTK_QUAD, 4, pts);
144 cell_code->SetValue(id_cell, part_bc[face[0]]);
145 orgdir->SetValue(id_cell, 0);
146 voldir->SetValue(id_cell, 0);
147 curdir->SetValue(id_cell, 0);
150 UpdateNodeIndex(m_Grid);
151 UpdateCellIndex(m_Grid);
153 // set the boundary names
154 GuiMainWindow::pointer()->clearBCs();
155 for (int i_part = 0; i_part < part_name.size(); ++i_part) {
156 GuiMainWindow::pointer()->addBC(part_bc[i_part], BoundaryCondition(part_name[i_part], "patch"));
160 } catch (Error err) {
161 err.display();