limited volume meshing to boundary layer only
[engrid-github.git] / src / libengrid / vtkEgBoundaryCodesFilter.cxx
blob58fb3083f0c6bad8dde6716cba35bd9611bd7ef0
1 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 // + +
3 // + This file is part of enGrid. +
4 // + +
5 // + Copyright 2008-2014 enGits GmbH +
6 // + +
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. +
11 // + +
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. +
16 // + +
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/>. +
19 // + +
20 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
21 #include "vtkEgBoundaryCodesFilter.h"
22 #include "engrid.h"
24 #include <vtkObjectFactory.h>
25 #include <vtkInformation.h>
26 #include <vtkDataObject.h>
27 #include <vtkPolyData.h>
28 #include <vtkUnstructuredGrid.h>
29 #include <vtkSmartPointer.h>
30 #include <vtkCellData.h>
31 #include <vtkIntArray.h>
33 vtkStandardNewMacro(vtkEgBoundaryCodesFilter)
35 void vtkEgBoundaryCodesFilter::ExecuteEg()
37 EG_VTKDCC(vtkIntArray,cell_code,m_Input,"cell_code");
39 // copy the node coordinates
40 allocateGrid(m_Output, m_Input->GetNumberOfCells(), m_Input->GetNumberOfPoints());
41 for (vtkIdType vertId = 0; vertId < m_Input->GetNumberOfPoints(); ++vertId) {
42 double x[3];
43 m_Input->GetPoints()->GetPoint(vertId,x);
44 m_Output->GetPoints()->SetPoint(vertId,x);
47 // copy the cells and the cell/node data
48 for (vtkIdType cellId = 0; cellId < m_Input->GetNumberOfCells(); ++cellId) {
49 bool add = false;
50 if (!cell_code) {
51 add = false;
53 else {
54 if (m_BoundaryCodes.contains(cell_code->GetValue(cellId))) {
55 if (m_Input->GetCellType(cellId) == VTK_TRIANGLE) add = true;
56 if (m_Input->GetCellType(cellId) == VTK_QUAD) add = true;
57 if (m_Input->GetCellType(cellId) == VTK_POLYGON) add = true;
61 if (add) {
62 EG_GET_CELL(cellId, m_Input);
63 vtkIdType newCell = m_Output->InsertNextCell(m_Input->GetCellType(cellId), ptIds);
64 copyCellData(m_Input, cellId, m_Output, newCell);
65 for (int i = 0; i < num_pts; i++) {
66 if (pts[i] >= m_Input->GetNumberOfPoints()) {
67 EG_BUG;
68 } else if (pts[i] >= m_Output->GetNumberOfPoints()) {
69 EG_BUG;
70 } else {
71 copyNodeData(m_Input, pts[i], m_Output, pts[i]);
77 m_Output->Squeeze();