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