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