Salome HOME
Merge with OCC_development_01
[modules/visu.git] / src / CONVERTOR / VISU_ExtractUnstructuredGrid.cxx
1 //  VISU CONVERTOR :
2 //
3 //  Copyright (C) 2003  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 //  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS 
5 // 
6 //  This library is free software; you can redistribute it and/or 
7 //  modify it under the terms of the GNU Lesser General Public 
8 //  License as published by the Free Software Foundation; either 
9 //  version 2.1 of the License. 
10 // 
11 //  This library is distributed in the hope that it will be useful, 
12 //  but WITHOUT ANY WARRANTY; without even the implied warranty of 
13 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
14 //  Lesser General Public License for more details. 
15 // 
16 //  You should have received a copy of the GNU Lesser General Public 
17 //  License along with this library; if not, write to the Free Software 
18 //  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA 
19 // 
20 //  See http://www.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org 
21 //
22 //
23 // File:    VISU_ExtractUnstructuredGrid.cxx
24 // Author:  Alexey PETROV
25 // Module : VISU
26
27
28 #include "VISU_ExtractUnstructuredGrid.hxx"
29 #include "VISU_ConvertorUtils.hxx"
30
31 #include <vtkUnstructuredGrid.h>
32 #include <vtkObjectFactory.h>
33 #include <vtkIdList.h>
34 #include <vtkCell.h>
35
36 using namespace std;
37
38 #ifdef _DEBUG_
39 static int MYDEBUG = 0;
40 #else
41 static int MYDEBUG = 0;
42 #endif
43
44 vtkStandardNewMacro(VISU_ExtractUnstructuredGrid);
45
46 VISU_ExtractUnstructuredGrid::VISU_ExtractUnstructuredGrid(){}
47
48 VISU_ExtractUnstructuredGrid::~VISU_ExtractUnstructuredGrid(){}
49
50 void VISU_ExtractUnstructuredGrid::RemoveCell(vtkIdType theCellId){
51   myRemovedCellIds.insert(theCellId);
52   Modified();
53 }
54
55 void VISU_ExtractUnstructuredGrid::RemoveCellsWithType(vtkIdType theCellType){
56   myRemovedCellTypes.insert(theCellType);
57   Modified();
58 }
59
60 namespace{
61   inline void InsertCell(vtkUnstructuredGrid *theInput, 
62                           vtkUnstructuredGrid *theOutput,
63                           vtkIdType theCellId, vtkIdList *theCellIds)
64   {
65     theCellIds->Reset();
66     vtkCell *aCell = theInput->GetCell(theCellId);
67     vtkIdType aNbIds = aCell->PointIds->GetNumberOfIds();
68     for(vtkIdType i = 0; i < aNbIds; i++)
69       theCellIds->InsertNextId(aCell->GetPointIds()->GetId(i));
70     theOutput->InsertNextCell(theInput->GetCellType(theCellId), theCellIds);
71   }
72 }
73
74 void VISU_ExtractUnstructuredGrid::Execute(){
75   vtkUnstructuredGrid *anInput = this->GetInput(), *anOutput = this->GetOutput();
76   vtkIdType aNbCells = anInput->GetNumberOfCells();
77   anOutput->Allocate(aNbCells);
78   MSG(MYDEBUG,"Execute - anInput->GetNumberOfCells() = "<<anInput->GetNumberOfCells());
79   vtkIdList *aCellIds = vtkIdList::New();
80   MSG(MYDEBUG,"Execute - myRemovedCellIds.empty() = "<<myRemovedCellIds.empty()<<
81       "; myRemovedCellTypes.empty() = "<<myRemovedCellTypes.empty());
82   if(myRemovedCellIds.empty() && myRemovedCellTypes.empty())
83     anOutput->CopyStructure(anInput);
84   else if(!myRemovedCellIds.empty() && myRemovedCellTypes.empty()){
85     for(vtkIdType aCellId = 0; aCellId < aNbCells; aCellId++)
86       if(myRemovedCellIds.find(aCellId) == myRemovedCellIds.end())
87         InsertCell(anInput,anOutput,aCellId,aCellIds);
88   }else if(myRemovedCellIds.empty() && !myRemovedCellTypes.empty()){
89     for(vtkIdType aCellId = 0; aCellId < aNbCells; aCellId++)
90       if(myRemovedCellTypes.find(anInput->GetCellType(aCellId)) == myRemovedCellTypes.end())
91         InsertCell(anInput,anOutput,aCellId,aCellIds);
92   }else if(!myRemovedCellIds.empty() && !myRemovedCellTypes.empty())
93     for(vtkIdType aCellId = 0; aCellId < aNbCells; aCellId++)
94       if(myRemovedCellTypes.find(anInput->GetCellType(aCellId)) == myRemovedCellTypes.end())
95         if(myRemovedCellIds.find(aCellId) == myRemovedCellIds.end())
96           InsertCell(anInput,anOutput,aCellId,aCellIds);
97   aCellIds->Delete();
98   anOutput->SetPoints(anInput->GetPoints());
99   MSG(MYDEBUG,"Execute - anOutput->GetNumberOfCells() = "<<anOutput->GetNumberOfCells());
100 }