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