Salome HOME
2c4293d6b10a00ed75ea1919d8a2bd1bd329941e
[modules/visu.git] / src / PIPELINE / SALOME_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 "SALOME_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
33 #include <vtkInformation.h>
34 #include <vtkInformationVector.h>
35
36 #include <vtkImplicitBoolean.h>
37 #include <vtkImplicitFunctionCollection.h>
38
39 #if defined __GNUC__
40   #if __GNUC__ == 2
41     #define __GNUC_2__
42   #endif
43 #endif
44
45
46 //----------------------------------------------------------------------------
47 vtkStandardNewMacro(SALOME_ExtractGeometry);
48
49
50 //----------------------------------------------------------------------------
51 SALOME_ExtractGeometry
52 ::SALOME_ExtractGeometry():
53   myStoreMapping(false),
54   myIsDoneShallowCopy(false)
55 {}
56
57 SALOME_ExtractGeometry
58 ::~SALOME_ExtractGeometry()
59 {}
60
61
62 //----------------------------------------------------------------------------
63 vtkImplicitBoolean* 
64 SALOME_ExtractGeometry
65 ::GetImplicitBoolean() 
66 {
67   return myImplicitBoolean.GetPointer();
68 }
69
70
71 void
72 SALOME_ExtractGeometry
73 ::SetImplicitFunction(vtkImplicitFunction* theImplicitFunction)  
74 {
75   myImplicitBoolean = dynamic_cast<vtkImplicitBoolean*>(theImplicitFunction);
76   Superclass::SetImplicitFunction(theImplicitFunction);
77 }
78
79
80 //----------------------------------------------------------------------------
81 void 
82 SALOME_ExtractGeometry
83 ::SetStoreMapping(bool theStoreMapping)
84 {
85   if(myStoreMapping != theStoreMapping){
86     myStoreMapping = theStoreMapping;
87     Modified();
88   }
89 }
90
91 bool 
92 SALOME_ExtractGeometry
93 ::GetStoreMapping() const
94 {
95   return myStoreMapping;
96 }
97
98
99 //----------------------------------------------------------------------------
100 vtkIdType
101 SALOME_ExtractGeometry
102 ::GetElemVTKId(vtkIdType theID)
103 {
104   if(!myStoreMapping || myIsDoneShallowCopy)
105     return theID;
106
107   vtkIdType iEnd = myElemVTK2ObjIds.size();
108   for(vtkIdType i = 0; i < iEnd; i++)
109     if(myElemVTK2ObjIds[i] == theID)
110       return i;
111
112   return -1;
113 }
114
115 vtkIdType
116 SALOME_ExtractGeometry
117 ::GetNodeVTKId(vtkIdType theID)
118 {
119   if(!myStoreMapping || myIsDoneShallowCopy)
120     return theID;
121
122   vtkIdType iEnd = myNodeVTK2ObjIds.size();
123   for(vtkIdType i = 0; i < iEnd; i++)
124     if(myNodeVTK2ObjIds[i] == theID)
125       return i;
126
127   return -1;
128 }
129
130
131 //----------------------------------------------------------------------------
132 vtkIdType
133 SALOME_ExtractGeometry
134 ::GetElemObjId(vtkIdType theVtkID)
135 {
136   if(!myStoreMapping || myIsDoneShallowCopy)
137     return theVtkID;
138
139   if(theVtkID < myElemVTK2ObjIds.size())
140     return myElemVTK2ObjIds[theVtkID];
141
142   return -1;
143 }
144
145
146 vtkIdType
147 SALOME_ExtractGeometry
148 ::GetNodeObjId(vtkIdType theVtkID)
149 {
150   if(!myStoreMapping || myIsDoneShallowCopy)
151     return theVtkID;
152
153   if(theVtkID < myNodeVTK2ObjIds.size())
154     return myNodeVTK2ObjIds[theVtkID];
155
156   return -1;
157 }
158
159
160 //----------------------------------------------------------------------------
161 int
162 SALOME_ExtractGeometry
163 ::RequestData(vtkInformation *request,
164               vtkInformationVector **inputVector,
165               vtkInformationVector *outputVector)
166 {
167   // get the info objects
168   vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
169   vtkInformation *outInfo = outputVector->GetInformationObject(0);
170
171   // get the input and ouptut
172   vtkDataSet *input = vtkDataSet::SafeDownCast(
173     inInfo->Get(vtkDataObject::DATA_OBJECT()));
174   vtkUnstructuredGrid *output = vtkUnstructuredGrid::SafeDownCast(
175     outInfo->Get(vtkDataObject::DATA_OBJECT()));
176
177   myElemVTK2ObjIds.clear();
178   myNodeVTK2ObjIds.clear();
179   //
180   myIsDoneShallowCopy = !this->ImplicitFunction;
181
182   if(!myIsDoneShallowCopy && myImplicitBoolean.GetPointer())
183     if(vtkImplicitFunctionCollection* aFunction = myImplicitBoolean->GetFunction())
184       myIsDoneShallowCopy = aFunction->GetNumberOfItems() == 0;
185
186   if(myIsDoneShallowCopy){
187     output->ShallowCopy(input);
188     return 1;
189   }
190   
191   return RequestData2(request,inputVector,outputVector);
192 }
193
194
195 //----------------------------------------------------------------------------
196 int
197 SALOME_ExtractGeometry
198 ::RequestData2(vtkInformation *vtkNotUsed(request),
199                vtkInformationVector **inputVector,
200                vtkInformationVector *outputVector)
201 {
202   // get the info objects
203   vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
204   vtkInformation *outInfo = outputVector->GetInformationObject(0);
205
206   // get the input and ouptut
207   vtkDataSet *input = vtkDataSet::SafeDownCast(
208     inInfo->Get(vtkDataObject::DATA_OBJECT()));
209   vtkUnstructuredGrid *output = vtkUnstructuredGrid::SafeDownCast(
210     outInfo->Get(vtkDataObject::DATA_OBJECT()));
211
212   vtkIdType ptId, numPts, numCells, i, cellId, newCellId, newId, *pointMap;
213   vtkIdList *cellPts;
214   vtkCell *cell;
215   int numCellPts;
216   vtkFloatingPointType *x;
217   vtkFloatingPointType multiplier;
218   vtkPoints *newPts;
219   vtkIdList *newCellPts;
220   vtkPointData *pd = input->GetPointData();
221   vtkCellData *cd = input->GetCellData();
222   vtkPointData *outputPD = output->GetPointData();
223   vtkCellData *outputCD = output->GetCellData();
224   int npts;
225   numCells = input->GetNumberOfCells();
226   numPts = input->GetNumberOfPoints();
227
228   if ( ! this->ImplicitFunction )
229     {
230     vtkErrorMacro(<<"No implicit function specified");
231     return 0;
232     }
233
234   newCellPts = vtkIdList::New();
235   newCellPts->Allocate(VTK_CELL_SIZE);
236
237   if ( this->ExtractInside )
238     {
239     multiplier = 1.0;
240     }
241   else
242     {
243     multiplier = -1.0;
244     }
245
246   // Loop over all points determining whether they are inside the
247   // implicit function. Copy the points and point data if they are.
248   //
249   pointMap = new vtkIdType[numPts]; // maps old point ids into new
250   for (i=0; i < numPts; i++)
251     {
252     pointMap[i] = -1;
253     }
254
255   output->Allocate(numCells/4); //allocate storage for geometry/topology
256   newPts = vtkPoints::New();
257   newPts->Allocate(numPts/4,numPts);
258   outputPD->CopyAllocate(pd);
259   outputCD->CopyAllocate(cd);
260   vtkFloatArray *newScalars = NULL;
261
262   if(myStoreMapping){
263     myElemVTK2ObjIds.reserve(numCells);
264     myNodeVTK2ObjIds.reserve(numPts);
265   }
266
267   if ( ! this->ExtractBoundaryCells )
268     {
269     for ( ptId=0; ptId < numPts; ptId++ )
270       {
271       x = input->GetPoint(ptId);
272       if ( (this->ImplicitFunction->FunctionValue(x)*multiplier) < 0.0 )
273         {
274         newId = newPts->InsertNextPoint(x);
275         pointMap[ptId] = newId;
276         if(myStoreMapping)
277           myNodeVTK2ObjIds.push_back(ptId);
278         outputPD->CopyData(pd,ptId,newId);
279         }
280       }
281     }
282   else
283     {
284     // To extract boundary cells, we have to create supplemental information
285     if ( this->ExtractBoundaryCells )
286       {
287       vtkFloatingPointType val;
288       newScalars = vtkFloatArray::New();
289       newScalars->SetNumberOfValues(numPts);
290
291       for (ptId=0; ptId < numPts; ptId++ )
292         {
293         x = input->GetPoint(ptId);
294         val = this->ImplicitFunction->FunctionValue(x) * multiplier;
295         newScalars->SetValue(ptId, val);
296         if ( val < 0.0 )
297           {
298           newId = newPts->InsertNextPoint(x);
299           pointMap[ptId] = newId;
300           if(myStoreMapping)
301             myNodeVTK2ObjIds.push_back(ptId);
302           outputPD->CopyData(pd,ptId,newId);
303           }
304         }
305       }
306     }
307
308   // Now loop over all cells to see whether they are inside implicit
309   // function (or on boundary if ExtractBoundaryCells is on).
310   //
311   for (cellId=0; cellId < numCells; cellId++)
312     {
313     cell = input->GetCell(cellId);
314     cellPts = cell->GetPointIds();
315     numCellPts = cell->GetNumberOfPoints();
316
317     newCellPts->Reset();
318     if ( ! this->ExtractBoundaryCells ) //requires less work
319       {
320       for ( npts=0, i=0; i < numCellPts; i++, npts++)
321         {
322         ptId = cellPts->GetId(i);
323         if ( pointMap[ptId] < 0 )
324           {
325           break; //this cell won't be inserted
326           }
327         else
328           {
329           newCellPts->InsertId(i,pointMap[ptId]);
330           }
331         }
332       } //if don't want to extract boundary cells
333
334     else //want boundary cells
335       {
336       for ( npts=0, i=0; i < numCellPts; i++ )
337         {
338         ptId = cellPts->GetId(i);
339         if ( newScalars->GetValue(ptId) <= 0.0 )
340           {
341           npts++;
342           }
343         }
344       if ( npts > 0 )
345         {
346         for ( i=0; i < numCellPts; i++ )
347           {
348           ptId = cellPts->GetId(i);
349           if ( pointMap[ptId] < 0 )
350             {
351             x = input->GetPoint(ptId);
352             newId = newPts->InsertNextPoint(x);
353             pointMap[ptId] = newId;
354             if(myStoreMapping)
355               myNodeVTK2ObjIds.push_back(ptId);
356             outputPD->CopyData(pd,ptId,newId);
357             }
358           newCellPts->InsertId(i,pointMap[ptId]);
359           }
360         }//a boundary or interior cell
361       }//if mapping boundary cells
362
363     if ( npts >= numCellPts || (this->ExtractBoundaryCells && npts > 0) )
364       {
365       newCellId = output->InsertNextCell(cell->GetCellType(),newCellPts);
366       if(myStoreMapping)
367         myElemVTK2ObjIds.push_back(cellId);
368       outputCD->CopyData(cd,cellId,newCellId);
369       }
370     }//for all cells
371
372   // Update ourselves and release memory
373   //
374   delete [] pointMap;
375   newCellPts->Delete();
376   output->SetPoints(newPts);
377   newPts->Delete();
378
379   if ( this->ExtractBoundaryCells )
380     {
381     newScalars->Delete();
382     }
383
384   output->Squeeze();
385
386   return 1;
387 }