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