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