2 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
4 // + This file is part of enGrid. +
6 // + Copyright 2008 Oliver Gloth +
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. +
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. +
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/>. +
21 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
23 #include "vtkEgGridFilter.h"
24 #include "vtkInformation.h"
26 #include "guimainwindow.h"
28 vtkEgGridFilter::vtkEgGridFilter()
30 m_BoundaryCodes
.clear();
33 m_LastOutput
= vtkUnstructuredGrid::New();
34 m_LastRun
= GetMTime();
37 vtkEgGridFilter::~vtkEgGridFilter()
39 m_LastOutput
->Delete();
42 int vtkEgGridFilter::FillInputPortInformation
48 info
->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(),"vtkUnstructuredGrid");
52 int vtkEgGridFilter::FillOutputPortInformation
58 info
->Set(vtkDataObject::DATA_TYPE_NAME(),"vtkUnstructuredGrid");
62 int vtkEgGridFilter::RequestData
65 vtkInformationVector
** m_InputVector
,
66 vtkInformationVector
* m_OutputVector
69 // Get m_Input and m_Output data
70 m_Input
= vtkUnstructuredGrid::GetData(m_InputVector
[0]);
71 m_Output
= vtkUnstructuredGrid::GetData(m_OutputVector
);
73 // check if we have a proper m_InputVector
74 if (!m_Input
) { vtkDebugMacro("No input!"); return 1; }
75 if (!m_Input
->GetPoints()) { vtkDebugMacro("No input!"); return 1; }
76 if (m_Input
->GetPoints()->GetNumberOfPoints() < 1) { vtkDebugMacro("No input!"); return 1; }
78 // call the enGrid filter method
80 bool m_Input_changed
= (m_Input
->GetMTime() > m_LastRun
);
81 bool filter_changed
= (GetMTime() > m_LastRun
);
82 if (m_Input_changed
|| filter_changed
) {
83 if (GuiMainWindow::tryLock()) {
85 makeCopy(m_Output
, m_LastOutput
);
87 m_LastRun
= GetMTime();
88 GuiMainWindow::unlock();
90 makeCopy(m_LastOutput
, m_Output
);
94 GuiMainWindow::unlock();
96 m_Output
->DeepCopy(m_Input
);
102 void vtkEgGridFilter::SetBoundaryCodes(const QSet
<int> &bc
)
105 if (m_BoundaryCodes
.size() != bc
.size()) {
108 QSet
<int> bc_inters
= m_BoundaryCodes
.intersect(bc
);
109 if (bc_inters
.size() != bc
.size()) {
114 m_BoundaryCodes
= bc
;
119 void vtkEgGridFilter::ExtractBoundary(QVector
<vtkIdType
> &cells
, QVector
<vtkIdType
> &nodes
, QSet
<int> &bc
, vtkUnstructuredGrid
*grid
)
123 QSet
<vtkIdType
> ex_nodes
, ex_cells
;
124 EG_VTKDCC(vtkIntArray
, cell_code
, grid
, "cell_code");
125 for (vtkIdType i
= 0; i
< grid
->GetNumberOfCells(); ++i
) {
126 if ((grid
->GetCellType(i
) == VTK_TRIANGLE
) || (grid
->GetCellType(i
) == VTK_QUAD
)){
127 if (bc
.contains(cell_code
->GetValue(i
))) {
131 m_Input
->GetCellPoints(i
,npts
,pts
);
132 for (int j
= 0; j
< npts
; ++j
) {
133 ex_nodes
.insert(pts
[j
]);
138 cells
.resize(ex_cells
.size());
139 nodes
.resize(ex_nodes
.size());
143 foreach(i
,ex_cells
) {
151 foreach(i
,ex_nodes
) {