Salome HOME
Merge from V5_1_3_BR branch (07/12/09)
[modules/smesh.git] / src / OBJECT / SMESH_ExtractGeometry.cxx
1 //  Copyright (C) 2007-2008  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 //  Copyright (C) 2003-2007  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.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 #include "SMESH_ExtractGeometry.h"
23
24 #include <vtkCell.h>
25 #include <vtkCellData.h>
26 #include <vtkFloatArray.h>
27 #include <vtkIdList.h>
28 #include <vtkImplicitFunction.h>
29 #include <vtkObjectFactory.h>
30 #include <vtkPointData.h>
31 #include <vtkUnstructuredGrid.h>
32 #include <vtkInformation.h>
33 #include <vtkInformationVector.h>
34
35 using namespace std;
36
37 #ifdef _DEBUG_
38 static int MYDEBUG = 0;
39 #else
40 static int MYDEBUG = 0;
41 #endif
42
43 #if defined __GNUC__
44   #if __GNUC__ == 2
45     #define __GNUC_2__
46   #endif
47 #endif
48
49
50 vtkStandardNewMacro(SMESH_ExtractGeometry);
51
52
53 SMESH_ExtractGeometry::SMESH_ExtractGeometry()
54 {}
55
56
57 SMESH_ExtractGeometry::~SMESH_ExtractGeometry(){}
58
59
60 vtkIdType SMESH_ExtractGeometry::GetElemObjId(int theVtkID){
61   if( theVtkID < 0 || theVtkID >= myElemVTK2ObjIds.size()) return -1;
62   return myElemVTK2ObjIds[theVtkID];
63 }
64
65
66 vtkIdType SMESH_ExtractGeometry::GetNodeObjId(int theVtkID){
67   if ( theVtkID < 0 || theVtkID >= myNodeVTK2ObjIds.size()) return -1;
68   return myNodeVTK2ObjIds[theVtkID];
69 }
70
71
72 int SMESH_ExtractGeometry::RequestData(
73   vtkInformation *vtkNotUsed(request),
74   vtkInformationVector **inputVector,
75   vtkInformationVector *outputVector)
76 {
77   // get the info objects
78   vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
79   vtkInformation *outInfo = outputVector->GetInformationObject(0);
80
81   // get the input and ouptut
82   vtkDataSet *input = vtkDataSet::SafeDownCast(
83     inInfo->Get(vtkDataObject::DATA_OBJECT()));
84   vtkUnstructuredGrid *output = vtkUnstructuredGrid::SafeDownCast(
85     outInfo->Get(vtkDataObject::DATA_OBJECT()));
86
87   vtkIdType ptId, numPts, numCells, i, cellId, newCellId, newId, *pointMap;
88   vtkIdList *cellPts;
89   vtkCell *cell;
90   int numCellPts;
91   vtkFloatingPointType *x;
92   vtkFloatingPointType multiplier;
93   vtkPoints *newPts;
94   vtkIdList *newCellPts;
95   vtkPointData *pd = input->GetPointData();
96   vtkCellData *cd = input->GetCellData();
97   vtkPointData *outputPD = output->GetPointData();
98   vtkCellData *outputCD = output->GetCellData();
99   int npts;
100   numCells = input->GetNumberOfCells();
101   numPts = input->GetNumberOfPoints();
102   
103   vtkDebugMacro(<< "Extracting geometry");
104
105   if ( ! this->ImplicitFunction )
106     {
107     vtkErrorMacro(<<"No implicit function specified");
108     return 0;
109     }
110
111   newCellPts = vtkIdList::New();
112   newCellPts->Allocate(VTK_CELL_SIZE);
113
114   if ( this->ExtractInside )
115     {
116     multiplier = 1.0;
117     }
118   else 
119     {
120     multiplier = -1.0;
121     }
122
123   // Loop over all points determining whether they are inside the
124   // implicit function. Copy the points and point data if they are.
125   //
126   pointMap = new vtkIdType[numPts]; // maps old point ids into new
127   for (i=0; i < numPts; i++)
128     {
129     pointMap[i] = -1;
130     }
131
132   output->Allocate(numCells/4); //allocate storage for geometry/topology
133   newPts = vtkPoints::New();
134   newPts->Allocate(numPts/4,numPts);
135   outputPD->CopyAllocate(pd);
136   outputCD->CopyAllocate(cd);
137   vtkFloatArray *newScalars = NULL;
138   
139   if(myStoreMapping){
140     myElemVTK2ObjIds.clear();
141     myElemVTK2ObjIds.reserve(numCells);
142     myNodeVTK2ObjIds.clear();
143     myNodeVTK2ObjIds.reserve(numPts);
144   }
145
146   if ( ! this->ExtractBoundaryCells )
147     {
148     for ( ptId=0; ptId < numPts; ptId++ )
149       {
150       x = input->GetPoint(ptId);
151       if ( (this->ImplicitFunction->FunctionValue(x)*multiplier) < 0.0 )
152         {
153         newId = newPts->InsertNextPoint(x);
154         pointMap[ptId] = newId;
155         myNodeVTK2ObjIds.push_back(ptId);
156         outputPD->CopyData(pd,ptId,newId);
157         }
158       }
159     }
160   else
161     {
162     // To extract boundary cells, we have to create supplemental information
163     if ( this->ExtractBoundaryCells )
164       {
165       vtkFloatingPointType val;
166       newScalars = vtkFloatArray::New();
167       newScalars->SetNumberOfValues(numPts);
168
169       for (ptId=0; ptId < numPts; ptId++ )
170         {
171         x = input->GetPoint(ptId);
172         val = this->ImplicitFunction->FunctionValue(x) * multiplier;
173         newScalars->SetValue(ptId, val);
174         if ( val < 0.0 )
175           {
176           newId = newPts->InsertNextPoint(x);
177           pointMap[ptId] = newId;
178           myNodeVTK2ObjIds.push_back(ptId);
179           outputPD->CopyData(pd,ptId,newId);
180           }
181         }
182       }
183     }
184
185   // Now loop over all cells to see whether they are inside implicit
186   // function (or on boundary if ExtractBoundaryCells is on).
187   //
188   for (cellId=0; cellId < numCells; cellId++)
189     {
190     cell = input->GetCell(cellId);
191     cellPts = cell->GetPointIds();
192     numCellPts = cell->GetNumberOfPoints();
193
194     newCellPts->Reset();
195     if ( ! this->ExtractBoundaryCells ) //requires less work
196       {
197       for ( npts=0, i=0; i < numCellPts; i++, npts++)
198         {
199         ptId = cellPts->GetId(i);
200         if ( pointMap[ptId] < 0 )
201           {
202           break; //this cell won't be inserted
203           }
204         else
205           {
206           newCellPts->InsertId(i,pointMap[ptId]);
207           }
208         }
209       } //if don't want to extract boundary cells
210     
211     else //want boundary cells
212       {
213       for ( npts=0, i=0; i < numCellPts; i++ )
214         {
215         ptId = cellPts->GetId(i);
216         if ( newScalars->GetValue(ptId) <= 0.0 )
217           {
218           npts++;
219           }
220         }
221       if ( npts > 0 )
222         {
223         for ( i=0; i < numCellPts; i++ )
224           {
225           ptId = cellPts->GetId(i);
226           if ( pointMap[ptId] < 0 )
227             {
228             x = input->GetPoint(ptId);
229             newId = newPts->InsertNextPoint(x);
230             pointMap[ptId] = newId;
231             myNodeVTK2ObjIds.push_back(ptId);
232             outputPD->CopyData(pd,ptId,newId);
233             }
234           newCellPts->InsertId(i,pointMap[ptId]);
235           }
236         }//a boundary or interior cell
237       }//if mapping boundary cells
238       
239     if ( npts >= numCellPts || (this->ExtractBoundaryCells && npts > 0) )
240       {
241       newCellId = output->InsertNextCell(cell->GetCellType(),newCellPts);
242       myElemVTK2ObjIds.push_back(cellId);
243       outputCD->CopyData(cd,cellId,newCellId);
244       }
245     }//for all cells
246
247   // Update ourselves and release memory
248   //
249   delete [] pointMap;
250   newCellPts->Delete();
251   output->SetPoints(newPts);
252   newPts->Delete();
253   
254   if ( this->ExtractBoundaryCells )
255     {
256     newScalars->Delete();
257     }
258
259   output->Squeeze();
260   return 1;
261 }