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 "setboundarycode.h"
22 #include <vtkIntArray.h>
23 #include <vtkCellData.h>
24 #include "guimainwindow.h"
26 SetBoundaryCode::SetBoundaryCode()
28 feature_angle
= 180.0;
30 setSurfaceIteration();
34 void SetBoundaryCode::pass1()
36 if(!(SelectAllVisible
|| OnlyPickedCell
|| OnlyPickedCellAndNeighbours
)) {
37 using namespace GeometryTools
;
38 double fa
= feature_angle
*M_PI
/180.0;
41 GuiMainWindow::pointer()->getDisplayBoundaryCodes(DBC
);
43 EG_VTKDCC(vtkIntArray
, cell_code
, m_Grid
, "cell_code");
45 for (int i
= 0; i
< pair
.size(); ++i
) {
46 int bc1
= cell_code
->GetValue(pair
[i
].item1
);
47 int bc2
= cell_code
->GetValue(pair
[i
].item2
);
49 vec3_t n1
= cellNormal(m_Grid
, pair
[i
].item1
);
50 vec3_t n2
= cellNormal(m_Grid
, pair
[i
].item2
);
53 if (GeometryTools::angle(n1
, n2
) > fa
) {
54 pair
[i
].terminate
= true;
56 pair
[i
].terminate
= false;
59 if (DBC
.contains(bc2
)) {
60 vec3_t n1
= cellNormal(m_Grid
, pair
[i
].item1
);
61 vec3_t n2
= cellNormal(m_Grid
, pair
[i
].item2
);
64 if (GeometryTools::angle(n1
, n2
) > fa
) {
65 pair
[i
].terminate
= true;
67 pair
[i
].terminate
= false;
70 pair
[i
].terminate
= true;
77 void SetBoundaryCode::pass2()
79 EG_VTKDCC(vtkIntArray
, cell_code
, m_Grid
, "cell_code");
81 if(SelectAllVisible
) {
83 GuiMainWindow::pointer()->getDisplayBoundaryCodes(DBC
);
84 DBC
.insert(boundary_code
);
85 foreach(cellId
, cells
) {
86 int bc
= cell_code
->GetValue(cellId
);
88 cell_code
->SetValue(cellId
, boundary_code
);
91 } else if (OnlyPickedCell
) {
92 cout
<<"this->getStart()="<<this->getStart()<<endl
;
93 cell_code
->SetValue(this->getStart(), boundary_code
);
94 } else if (OnlyPickedCellAndNeighbours
) {
95 cell_code
->SetValue(this->getStart(), boundary_code
);
96 foreach (cellId
, c2c
[this->getStart()]) {
97 cell_code
->SetValue(cellId
, boundary_code
);
100 foreach (cellId
, item
) {
101 cell_code
->SetValue(cellId
, boundary_code
);