]> SALOME platform Git repositories - modules/gui.git/blob - src/VTKViewer/VTKViewer_ExtractUnstructuredGrid.cxx
Salome HOME
23315: [CEA 1929] Too much memory used to display a mesh in shading and wireframe
[modules/gui.git] / src / VTKViewer / VTKViewer_ExtractUnstructuredGrid.cxx
1 // Copyright (C) 2007-2016  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 // File:    VTKViewer_ExtractUnstructuredGrid.cxx
24 // Author:  Alexey PETROV
25
26 #include "VTKViewer_ExtractUnstructuredGrid.h"
27 #include "VTKViewer_CellLocationsArray.h"
28
29 #include <vtkUnsignedCharArray.h>
30 #include <vtkUnstructuredGrid.h>
31 #include <vtkObjectFactory.h>
32 #include <vtkCellArray.h>
33 #include <vtkIdList.h>
34 #include <vtkCell.h>
35 #include <vtkCellData.h>
36 #include <vtkInformation.h>
37 #include <vtkInformationVector.h>
38 #include <vtkVersion.h>
39
40 #include "utilities.h"
41
42 #if defined __GNUC__
43   #if __GNUC__ == 2
44     #define __GNUC_2__
45   #endif
46 #endif
47
48 #define VTK_XVERSION (VTK_MAJOR_VERSION*10000+VTK_MINOR_VERSION*100+VTK_BUILD_VERSION)
49
50 vtkStandardNewMacro(VTKViewer_ExtractUnstructuredGrid);
51
52
53 VTKViewer_ExtractUnstructuredGrid::VTKViewer_ExtractUnstructuredGrid():
54   myExtractionMode(eCells), myChangeMode(ePassAll)
55 {}
56
57
58 VTKViewer_ExtractUnstructuredGrid::~VTKViewer_ExtractUnstructuredGrid()
59 {}
60
61 void VTKViewer_ExtractUnstructuredGrid::RegisterCell(vtkIdType theCellId)
62 {
63   if ( myCellIds.insert(theCellId).second )
64     Modified();
65 }
66
67
68 void VTKViewer_ExtractUnstructuredGrid::RegisterCellsWithType(vtkIdType theCellType)
69 {
70   if ( myCellTypes.insert(theCellType).second )
71     Modified();
72 }
73
74
75 void VTKViewer_ExtractUnstructuredGrid::SetStoreMapping(int theStoreMapping)
76 {
77   if ( myStoreMapping != ( theStoreMapping != 0 ))
78   {
79     myStoreMapping = theStoreMapping != 0;
80     Modified();
81   }
82 }
83
84 vtkIdType VTKViewer_ExtractUnstructuredGrid::GetInputId(int theOutId) const
85 {
86   if ( myCellIds.empty() && myCellTypes.empty() )
87     return theOutId;
88
89   if ( theOutId<0 || theOutId >= (int)myOut2InId.size() )
90     return -1;
91   return myOut2InId[theOutId];
92 }
93
94 // vtkIdType VTKViewer_ExtractUnstructuredGrid::GetOutputId(int theInId) const{
95 //   if(myCellIds.empty() && myCellTypes.empty()) return theInId;
96 //   TMapId::const_iterator anIter = myIn2OutId.find(theInId);
97 //   if(anIter == myIn2OutId.end()) return -1;
98 //   return anIter->second;
99 // }
100
101
102 inline int InsertCell(vtkUnstructuredGrid *theInput,
103                       vtkCellArray *theConnectivity,
104                       vtkUnsignedCharArray* theCellTypesArray,
105                       vtkIdTypeArray*& theFaces,
106                       vtkIdTypeArray*& theFaceLocations,
107                       vtkIdType theCellId,
108                       vtkIdList *theIdList,
109                       bool theStoreMapping,
110                       vtkIdType theOutId,
111                       VTKViewer_ExtractUnstructuredGrid::TVectorId& theOut2InId/*,
112                       VTKViewer_ExtractUnstructuredGrid::TMapId& theIn2OutId*/)
113 {
114   vtkCell      *aCell = theInput->GetCell(theCellId);
115   vtkIdType aCellType = aCell->GetCellType();
116   vtkIdType  aCellId = -1;
117 #if VTK_XVERSION > 50700
118   if (aCellType != VTK_POLYHEDRON)
119   {
120 #endif
121     aCellId = theConnectivity->InsertNextCell( aCell->GetPointIds() ); //theIdList);
122     if (theFaceLocations)
123       theFaceLocations->InsertNextValue(-1);
124 #if VTK_XVERSION > 50700
125   }
126   else
127   {
128     //MESSAGE("InsertCell type VTK_POLYHEDRON " << theStoreMapping);
129     if (!theFaces)
130     {
131       theFaces = vtkIdTypeArray::New();
132       theFaces->Allocate(theCellTypesArray->GetSize());
133       theFaceLocations = vtkIdTypeArray::New();
134       theFaceLocations->Allocate(theCellTypesArray->GetSize());
135       // FaceLocations must be padded until the current position
136       for (vtkIdType i = 0; i <= theCellTypesArray->GetMaxId(); i++)
137       {
138         theFaceLocations->InsertNextValue(-1);
139       }
140     }
141     // insert face location
142     theFaceLocations->InsertNextValue(theFaces->GetMaxId() + 1);
143
144     // insert cell connectivity and faces stream
145     vtkIdType nfaces = 0;
146     vtkIdType*  face = 0;
147     vtkIdType realnpts;
148     theInput->GetFaceStream(theCellId, nfaces, face);
149     vtkUnstructuredGrid::DecomposeAPolyhedronCell(
150                                                   nfaces, face, realnpts, theConnectivity, theFaces);
151   }
152 #endif
153
154   /*vtkIdType anID = */theCellTypesArray->InsertNextValue(aCellType);
155   if(theStoreMapping){
156     theOut2InId.push_back( theCellId );
157     //theIn2OutId.insert( theIn2OutId.end(), std::make_pair( theCellId, theOutId ));
158   }
159   return aCellId;
160 }
161
162 inline void InsertPointCell(vtkCellArray *theConnectivity, 
163                             vtkUnsignedCharArray* theCellTypesArray,
164                             vtkIdType theCellId,
165                             vtkIdList *theIdList,
166                             bool theStoreMapping,
167                             vtkIdType theOutId, 
168                             VTKViewer_ExtractUnstructuredGrid::TVectorId& theOut2InId/*,
169                             VTKViewer_ExtractUnstructuredGrid::TMapId& theIn2OutId*/)
170 {
171   theIdList->SetId(0,theCellId);
172   theConnectivity->InsertNextCell(theIdList);
173   theCellTypesArray->InsertNextValue(VTK_VERTEX);
174   if(theStoreMapping){
175     theOut2InId.push_back(theCellId);
176     //theIn2OutId.insert( theIn2OutId.end(), std::make_pair( theCellId, theOutId ));
177   }
178 }
179
180
181 int VTKViewer_ExtractUnstructuredGrid::RequestData(vtkInformation *vtkNotUsed(request),
182                                                    vtkInformationVector **inputVector,
183                                                    vtkInformationVector *outputVector)
184 {
185   // get the info objects
186   vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
187   vtkInformation *outInfo = outputVector->GetInformationObject(0);
188
189   // get the input and ouptut
190   vtkUnstructuredGrid *anInput = vtkUnstructuredGrid::SafeDownCast(
191     inInfo->Get(vtkDataObject::DATA_OBJECT()));
192   vtkUnstructuredGrid *anOutput = vtkUnstructuredGrid::SafeDownCast(
193     outInfo->Get(vtkDataObject::DATA_OBJECT()));
194
195   //vtkUnstructuredGrid *anInput = this->GetInput();
196   //vtkUnstructuredGrid *anOutput = this->GetOutput();
197
198   myOut2InId.clear();  //myIn2OutId.clear();
199
200   // use a vector of cellTypes to avoid searching in myCellTypes map
201   // for a better performance (IPAL53103)
202   TVectorId cellTypesVec( VTK_NUMBER_OF_CELL_TYPES, -1 );
203   for ( TSetId::iterator type = myCellTypes.begin(); type != myCellTypes.end(); ++type )
204   {
205     if ( *type >= (int)cellTypesVec.size() ) cellTypesVec.resize( *type+1, -1 );
206     if ( *type > 0 )
207       cellTypesVec[ *type ] = *type;
208   }
209
210 /*  if(MYDEBUG){
211     MESSAGE("Execute - anInput->GetNumberOfCells() = "<<anInput->GetNumberOfCells());
212     MESSAGE("Execute - myCellTypes.size() = "<<myCellTypes.size());
213     MESSAGE("Execute - myCellIds.size() = "<<myCellIds.size());
214     MESSAGE("Execute - myExtractionMode = "<<myExtractionMode);
215     MESSAGE("Execute - myChangeMode = "<<myChangeMode);
216   }*/
217   if(myExtractionMode == eCells){
218     if(myChangeMode == ePassAll || (myCellIds.empty() && myCellTypes.empty() && myChangeMode == eRemoving)){
219       if(vtkIdType aNbElems = anInput->GetNumberOfCells()){
220         if(myStoreMapping) myOut2InId.reserve(aNbElems);
221         anOutput->ShallowCopy(anInput);
222         if(myStoreMapping){
223           for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
224             myOut2InId.push_back(aCellId);
225             //myIn2OutId.insert( myIn2OutId.end(), std::make_pair( aCellId, anOutId ));
226           }
227         }
228       }
229     }else{
230       vtkIdList *anIdList = vtkIdList::New();
231       vtkCellArray *aConnectivity = vtkCellArray::New();
232       vtkIdType aNbElems = anInput->GetNumberOfCells();
233       aConnectivity->Allocate(2*aNbElems,0);
234       vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
235       aCellTypesArray->SetNumberOfComponents(1);
236       aCellTypesArray->Allocate(aNbElems*aCellTypesArray->GetNumberOfComponents());
237       anOutput->GetCellData()->CopyAllocate(anInput->GetCellData(),aNbElems,aNbElems/2);
238
239       vtkIdTypeArray *newFaces = 0;
240       vtkIdTypeArray *newFaceLocations = 0;
241
242       if(!myCellIds.empty() && myCellTypes.empty()){
243         if(myStoreMapping) myOut2InId.reserve(myCellIds.size());
244         if(myChangeMode == eAdding){
245           for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
246             if(myCellIds.find(aCellId) != myCellIds.end()){
247               vtkIdType newId = InsertCell(anInput,aConnectivity,aCellTypesArray,newFaces,newFaceLocations,aCellId,anIdList,
248                                            myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
249               anOutput->GetCellData()->CopyData(anInput->GetCellData(),aCellId,newId);
250             }
251           }
252         }else{
253           for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
254             if(myCellIds.find(aCellId) == myCellIds.end()){
255               vtkIdType newId = InsertCell(anInput,aConnectivity,aCellTypesArray,newFaces,newFaceLocations,aCellId,anIdList,
256                                            myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
257               anOutput->GetCellData()->CopyData(anInput->GetCellData(),aCellId,newId);
258             }
259           }
260         }
261       }else if(myCellIds.empty() && !myCellTypes.empty()){
262         if(myChangeMode == eAdding){
263           for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
264             vtkIdType aType = anInput->GetCellType(aCellId);
265             if ( cellTypesVec[ aType ] == aType ) {
266               vtkIdType newId = InsertCell(anInput,aConnectivity,aCellTypesArray,newFaces,newFaceLocations,aCellId,anIdList,
267                                            myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
268               anOutput->GetCellData()->CopyData(anInput->GetCellData(),aCellId,newId);
269             }
270           }
271         }else{
272           for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
273             vtkIdType aType = anInput->GetCellType(aCellId);
274             if ( cellTypesVec[ aType ] != aType ) {
275               vtkIdType newId = InsertCell(anInput,aConnectivity,aCellTypesArray,newFaces,newFaceLocations,aCellId,anIdList,
276                                            myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
277               anOutput->GetCellData()->CopyData(anInput->GetCellData(),aCellId,newId);
278             }
279           }
280         }
281       }else if(!myCellIds.empty() && !myCellTypes.empty()){
282         if(myChangeMode == eAdding){
283           for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
284             vtkIdType aType = anInput->GetCellType(aCellId);
285             if ( cellTypesVec[ aType ] == aType ) {
286               if(myCellIds.find(aCellId) != myCellIds.end()){
287                 vtkIdType newId = InsertCell(anInput,aConnectivity,aCellTypesArray,newFaces,newFaceLocations,aCellId,anIdList,
288                                              myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
289                 anOutput->GetCellData()->CopyData(anInput->GetCellData(),aCellId,newId);
290               }
291             }
292           }
293         }else{
294           for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
295             vtkIdType aType = anInput->GetCellType(aCellId);
296             if ( cellTypesVec[ aType ] != aType ) {
297               if(myCellIds.find(aCellId) == myCellIds.end()){
298                 vtkIdType newId = InsertCell(anInput,aConnectivity,aCellTypesArray,newFaces,newFaceLocations,aCellId,anIdList,
299                                              myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
300                 anOutput->GetCellData()->CopyData(anInput->GetCellData(),aCellId,newId);
301               }
302             }
303           }
304         }
305       }
306       if((aNbElems = aConnectivity->GetNumberOfCells())){
307         VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
308         aCellLocationsArray->SetNumberOfComponents(1);
309         aCellLocationsArray->SetNumberOfTuples(aNbElems);
310         aConnectivity->InitTraversal();
311         for(vtkIdType i = 0, *pts, npts; aConnectivity->GetNextCell(npts,pts); i++){
312           aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
313         }
314 #if VTK_XVERSION > 50700
315         anOutput->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity,newFaceLocations,newFaces);
316 #else
317         anOutput->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
318 #endif
319         anOutput->SetPoints(anInput->GetPoints());
320         aCellLocationsArray->Delete();
321       }
322       aCellTypesArray->Delete();
323       aConnectivity->Delete();
324       anIdList->Delete();
325       if ( newFaceLocations ) newFaceLocations->Delete();
326       if ( newFaces ) newFaces->Delete();
327     }
328   }else{
329     vtkIdList *anIdList = vtkIdList::New();
330     anIdList->SetNumberOfIds(1);
331     vtkCellArray *aConnectivity = vtkCellArray::New();
332     vtkIdType aNbElems = anInput->GetNumberOfPoints();
333     aConnectivity->Allocate(2*aNbElems,0);
334     vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
335     aCellTypesArray->SetNumberOfComponents(1);
336     aCellTypesArray->Allocate(aNbElems*aCellTypesArray->GetNumberOfComponents());
337     // additional condition has been added to treat a case described in IPAL21372
338     // note that it is significant only when myExtractionMode == ePoints
339     if(myChangeMode == ePassAll || (myCellIds.empty() && myCellTypes.empty() && myChangeMode == eRemoving) ||
340        !anInput->GetCellTypesArray()){
341       if(myStoreMapping) myOut2InId.reserve(aNbElems);
342       for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
343         InsertPointCell(aConnectivity,aCellTypesArray,aCellId,anIdList,
344                         myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
345       }
346     }else if(!myCellIds.empty() && myCellTypes.empty()){
347       if(myStoreMapping) myOut2InId.reserve(myCellIds.size());
348       if(myChangeMode == eAdding){
349         for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
350           if(myCellIds.find(aCellId) != myCellIds.end()){
351             InsertPointCell(aConnectivity,aCellTypesArray,aCellId,anIdList,
352                             myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
353           }
354         }
355       }else{
356         for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
357           if(myCellIds.find(aCellId) == myCellIds.end()){
358             InsertPointCell(aConnectivity,aCellTypesArray,aCellId,anIdList,
359                             myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
360           }
361         }
362       }
363     }else if(myCellIds.empty() && !myCellTypes.empty()){
364       if(myChangeMode == eAdding){
365         for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
366           vtkIdType aType = anInput->GetCellType(aCellId);
367           if ( cellTypesVec[ aType ] == aType ) {
368             InsertPointCell(aConnectivity,aCellTypesArray,aCellId,anIdList,
369                             myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
370           }
371         }
372       }else{
373         for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
374           vtkIdType aType = anInput->GetCellType(aCellId);
375           if ( cellTypesVec[ aType ] != aType ) {
376             InsertPointCell(aConnectivity,aCellTypesArray,aCellId,anIdList,
377                             myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
378           }
379         }
380       }
381     }else if(!myCellIds.empty() && !myCellTypes.empty()){
382       if(myChangeMode == eAdding){
383         for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
384           vtkIdType aType = anInput->GetCellType(aCellId);
385           if ( cellTypesVec[ aType ] == aType ) {
386             if(myCellIds.find(aCellId) != myCellIds.end()){
387               InsertPointCell(aConnectivity,aCellTypesArray,aCellId,anIdList,
388                               myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
389             }
390           }
391         }
392       }else{
393         for(vtkIdType aCellId = 0, anOutId = 0; aCellId < aNbElems; aCellId++,anOutId++){
394           vtkIdType aType = anInput->GetCellType(aCellId);
395           if ( cellTypesVec[ aType ] != aType ) {
396             if(myCellIds.find(aCellId) == myCellIds.end()){
397               InsertPointCell(aConnectivity,aCellTypesArray,aCellId,anIdList,
398                               myStoreMapping,anOutId,myOut2InId/*,myIn2OutId*/);
399             }
400           }
401         }
402       }
403     }
404     if((aNbElems = aConnectivity->GetNumberOfCells())){
405       VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
406       aCellLocationsArray->SetNumberOfComponents(1);
407       aCellLocationsArray->SetNumberOfTuples(aNbElems);
408       aConnectivity->InitTraversal();
409       for(vtkIdType i = 0, *pts, npts; aConnectivity->GetNextCell(npts,pts); i++){
410         aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
411       }
412 #if VTK_XVERSION > 50700
413       anOutput->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity,0, 0);
414 #else
415       anOutput->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
416 #endif
417       anOutput->SetPoints(anInput->GetPoints());
418       aCellLocationsArray->Delete();
419     }
420     aCellTypesArray->Delete();
421     aConnectivity->Delete();
422     anIdList->Delete();
423   }
424 /*  if(MYDEBUG){
425     MESSAGE("Execute - anOutput->GetNumberOfCells() = "<<anOutput->GetNumberOfCells());
426     if(myStoreMapping){
427       MESSAGE("Execute - myOut2InId.size() = "<<myOut2InId.size());
428       MESSAGE("Execute - myIn2OutId.size() = "<<myIn2OutId.size());
429     }
430   }*/
431   return 1;
432 }