Salome HOME
Compatibility CMake
[modules/visu.git] / src / CONVERTOR / VISU_AppendFilterUtilities.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 //  SALOME OBJECT : kernel of SALOME component
23 //  File   : VISU_GeometryFilter.cxx
24 //  Author : 
25 //  Module : SALOME
26 //  $Header$
27 //
28 #include "VISU_AppendFilterUtilities.hxx"
29 #include "VISU_ConvertorUtils.hxx"
30
31 #include <vtkCell.h>
32 #include <vtkCellData.h>
33 #include <vtkPointData.h>
34
35 #include <vtkDataSetCollection.h>
36 #include <vtkObjectFactory.h>
37
38 #include <vtkUnstructuredGrid.h>
39 #include <vtkPolyData.h>
40
41 #include <vtkInformationVector.h>
42 #include <vtkInformation.h>
43 #include <vtkExecutive.h>
44
45 #include <vtkPoints.h>
46 #include <vtkIntArray.h>
47
48 #include <algorithm>
49 #include <vector>
50 #include <map>
51
52 namespace
53 {
54   //---------------------------------------------------------------
55   typedef vtkIdType TCellId;
56   typedef vtkIdType TInputId;
57   typedef std::pair<TInputId, TCellId> TInputCellId;
58
59   typedef std::pair<vtkIdType, vtkIdType> TObjectId; 
60   typedef std::map<TObjectId, TInputCellId> TObject2InputIdMap;
61   
62
63   //---------------------------------------------------------------
64   void
65   DoMergingInputs(vtkCellData *theCellData, 
66                   TInputId theInputId,
67                   TObject2InputIdMap& theResult)
68   {
69     if(vtkDataArray *aDataArray = theCellData->GetArray("VISU_CELLS_MAPPER")){
70       if(vtkIntArray *anIntArray = dynamic_cast<vtkIntArray*>(aDataArray)){
71         int *aPointer = anIntArray->GetPointer(0);
72         int aNbCells = anIntArray->GetNumberOfTuples();
73         for(vtkIdType aCellId = 0; aCellId < aNbCells; aCellId++){
74           int aObjId = *aPointer++;
75           int anEntity = *aPointer++;
76           TObjectId anObjectId(aObjId, anEntity);
77           TObject2InputIdMap::iterator anIter = theResult.find(anObjectId);
78           if(anIter != theResult.end())
79             continue;
80           TInputCellId anInputCellId(theInputId, aCellId);
81           theResult.insert(anIter, TObject2InputIdMap::value_type(anObjectId, anInputCellId));
82         }
83       }
84     }
85   }
86
87
88   //---------------------------------------------------------------
89   struct TFillFieldList
90   {
91     vtkDataSetAttributes::FieldList myFieldList;
92     bool myIsFirstCellData;
93
94     TFillFieldList(vtkIdType theNbInputs):
95       myFieldList(theNbInputs),
96       myIsFirstCellData(true)
97     {}
98
99     void
100     operator()(TInputId theInputId, vtkDataSet* theDataSet)
101     {
102       vtkCellData *aCellData = theDataSet->GetCellData();
103       if(myIsFirstCellData){
104         myFieldList.InitializeFieldList(aCellData);
105         myIsFirstCellData = false;
106       }else{
107         myFieldList.IntersectFieldList(aCellData);
108       }
109     }
110     
111     virtual
112     vtkIdType
113     GetNbCells() const = 0;
114   };
115   
116
117   //---------------------------------------------------------------
118   struct TCellCounter: TFillFieldList
119   {
120     vtkIdType myNbCells;
121     
122     TCellCounter(vtkIdType theNbInputs):
123       TFillFieldList(theNbInputs),
124       myNbCells(0)
125     {}
126
127     void
128     operator()(TInputId theInputId, vtkDataSet* theDataSet)
129     {
130       TFillFieldList::operator()(theInputId, theDataSet);
131       myNbCells += theDataSet->GetNumberOfCells();
132     }
133
134     virtual
135     vtkIdType
136     GetNbCells() const
137     {
138       return myNbCells;
139     }
140   };
141
142
143   //---------------------------------------------------------------
144   struct TCellIdMerger: TFillFieldList
145   {
146     TObject2InputIdMap myObject2InputIdMap;
147
148     TCellIdMerger(vtkIdType theNbInputs):
149       TFillFieldList(theNbInputs)
150     {}
151
152     void
153     operator()(TInputId theInputId, vtkDataSet* theDataSet)
154     {
155       TFillFieldList::operator()(theInputId, theDataSet);
156       vtkCellData *aCellData = theDataSet->GetCellData();
157       DoMergingInputs(aCellData, theInputId, myObject2InputIdMap);
158     }
159
160     virtual
161     vtkIdType
162     GetNbCells() const
163     {
164       return myObject2InputIdMap.size();
165     }
166   };
167
168
169   //---------------------------------------------------------------
170   template<class TFunctor>
171   void
172   ForEachInput(vtkInformationVector **theInputVector, 
173                vtkIdType theNumberOfInputConnections,
174                TFunctor& theFunctor)
175   {
176     for(vtkIdType anInputId = 0; anInputId < theNumberOfInputConnections; anInputId++)
177       if(vtkDataSet *aDataSet = VISU::GetInput(theInputVector, anInputId))
178         if(aDataSet->GetNumberOfPoints() > 0 && aDataSet->GetNumberOfCells() > 0)
179           theFunctor(anInputId, aDataSet);
180   }
181
182
183   //---------------------------------------------------------------
184   template<class TDataSet>
185   bool
186   RequestData(vtkInformationVector **theInputVector,
187               vtkIdType theNumberOfInputConnections,
188               vtkInformationVector *theOutputVector,
189               vtkPointSet* theSharedPointSet,
190               bool theIsMergingInputs,
191               bool theIsMappingInputs)
192   {
193     if ( theNumberOfInputConnections == 1 ) {
194       // get the input and ouptut
195       vtkDataSet *anInput = VISU::GetInput( theInputVector, 0 );
196       vtkDataSet* anOutput = VISU::GetOutput( theOutputVector );
197
198       if ( anInput->GetDataObjectType() != anOutput->GetDataObjectType() )
199         return false;
200
201       // This has to be here because it initialized all field datas.
202       anOutput->CopyStructure( anInput );
203   
204       // Pass all. (data object's field data is passed by the
205       // superclass after this method)
206       anOutput->GetPointData()->PassData( anInput->GetPointData() );
207       anOutput->GetCellData()->PassData( anInput->GetCellData() );
208       
209       return true;
210     }
211
212     if ( theSharedPointSet ) {
213       vtkPoints* aPoints = theSharedPointSet->GetPoints();
214       if(aPoints->GetNumberOfPoints() < 1)
215         return true;
216   
217       TDataSet* anOutput = TDataSet::SafeDownCast(VISU::GetOutput(theOutputVector));
218       vtkIdType anNbInputs = theNumberOfInputConnections;
219       if ( theIsMergingInputs ) {
220         TCellIdMerger aFunctor(anNbInputs);
221         ForEachInput<TCellIdMerger>(theInputVector, anNbInputs, aFunctor);
222
223         vtkDataSetAttributes::FieldList& aFieldList = aFunctor.myFieldList;
224         TObject2InputIdMap& anObject2InputIdMap = aFunctor.myObject2InputIdMap;
225         vtkIdType aNbCells = aFunctor.GetNbCells();
226         if(aNbCells < 1)
227           return true;
228     
229         // Now can allocate memory
230         anOutput->Allocate(aNbCells); 
231         vtkCellData *anOutputCellData = anOutput->GetCellData();
232         anOutputCellData->CopyAllocate(aFieldList, aNbCells);
233       
234         // Append each input dataset together
235         // 1.points
236         anOutput->SetPoints(theSharedPointSet->GetPoints());
237         anOutput->GetPointData()->PassData(theSharedPointSet->GetPointData());
238       
239         // 2.cells
240         vtkIdList *anIdList = vtkIdList::New(); 
241         anIdList->Allocate(VTK_CELL_SIZE);
242         TObject2InputIdMap::const_iterator anIter = anObject2InputIdMap.begin();
243         TObject2InputIdMap::const_iterator anEndIter = anObject2InputIdMap.end();
244         for(; anIter != anEndIter; anIter++){
245           //TObjectId anObjectId = anIter->first;
246           const TInputCellId& anInputCellId = anIter->second;
247           TInputId anInputId = anInputCellId.first;
248           if(vtkDataSet *aDataSet = VISU::GetInput(theInputVector, anInputId)){
249             TCellId aCellId = anInputCellId.second;
250             aDataSet->GetCellPoints(aCellId, anIdList);
251             
252             vtkIdType aCellType = aDataSet->GetCellType(aCellId);
253             vtkIdType aNewCellId = anOutput->InsertNextCell(aCellType, anIdList);
254           
255             vtkCellData *aCellData = aDataSet->GetCellData();
256             anOutputCellData->CopyData(aFieldList, aCellData, anInputId, aCellId, aNewCellId);
257           }
258         }
259         anIdList->Delete();
260
261         if(theIsMappingInputs){
262           vtkIntArray *aDataArray = vtkIntArray::New();
263           aDataArray->SetName("VISU_INPUTS_MAPPER");
264           aDataArray->SetNumberOfComponents(2);
265           aDataArray->SetNumberOfTuples(aNbCells);
266
267           vtkIdType aTupleId = 0;
268           TObject2InputIdMap::const_iterator anIter = anObject2InputIdMap.begin();
269           TObject2InputIdMap::const_iterator anEndIter = anObject2InputIdMap.end();
270           for(vtkIdType aCellId = 0; anIter != anEndIter; anIter++, aCellId++){
271             const TInputCellId& anInputCellId = anIter->second;
272             TInputId anInputId = anInputCellId.first;
273             /*TCellId*/ aCellId = anInputCellId.second;
274             aDataArray->SetValue(aTupleId++, anInputId);
275             aDataArray->SetValue(aTupleId++, aCellId);
276           }
277
278           anOutputCellData->AddArray(aDataArray);
279           aDataArray->Delete();
280         }
281
282         return true;
283       }else{
284         TCellCounter aFunctor(anNbInputs);
285         ForEachInput<TCellCounter>(theInputVector, anNbInputs, aFunctor);
286         
287         vtkDataSetAttributes::FieldList& aFieldList = aFunctor.myFieldList;
288         vtkIdType aNbCells = aFunctor.GetNbCells();
289         if(aNbCells < 1)
290           return true;
291         
292         // Now can allocate memory
293         anOutput->Allocate(aNbCells); 
294         vtkCellData *anOutputCellData = anOutput->GetCellData();
295         anOutputCellData->CopyAllocate(aFieldList, aNbCells);
296         
297         // Append each input dataset together
298         // 1.points
299         anOutput->SetPoints(theSharedPointSet->GetPoints());
300         anOutput->GetPointData()->PassData(theSharedPointSet->GetPointData());
301         
302         // 2.cells
303         vtkIdList *anIdList = vtkIdList::New(); 
304         anIdList->Allocate(VTK_CELL_SIZE);
305         for(vtkIdType anInputId = 0; anInputId < anNbInputs; anInputId++){
306           if(vtkDataSet *aDataSet = VISU::GetInput(theInputVector, anInputId)){
307             vtkIdType aNbCells = aDataSet->GetNumberOfCells(); 
308             vtkCellData *aCellData = aDataSet->GetCellData();
309             // copy cell and cell data
310             for(vtkIdType aCellId = 0; aCellId < aNbCells; aCellId++){
311               aDataSet->GetCellPoints(aCellId, anIdList);
312
313               vtkIdType aCellType = aDataSet->GetCellType(aCellId);
314               vtkIdType aNewCellId = anOutput->InsertNextCell(aCellType, anIdList);
315
316               anOutputCellData->CopyData(aFieldList, aCellData, anInputId, aCellId, aNewCellId);
317             }
318           }
319         }
320         anIdList->Delete();
321
322         if(theIsMappingInputs){
323           vtkIntArray *aDataArray = vtkIntArray::New();
324           aDataArray->SetName("VISU_INPUTS_MAPPER");
325           aDataArray->SetNumberOfComponents(2);
326           aDataArray->SetNumberOfTuples(aNbCells);
327
328           vtkIdType aTupleId = 0;
329           for(vtkIdType anInputId = 0; anInputId < anNbInputs; anInputId++){
330             if(vtkDataSet *aDataSet = VISU::GetInput(theInputVector, anInputId)){
331               vtkIdType aNbCells = aDataSet->GetNumberOfCells(); 
332               for(vtkIdType aCellId = 0; aCellId < aNbCells; aCellId++){
333                 aDataArray->SetValue(aTupleId++, aCellId);
334                 aDataArray->SetValue(aTupleId++, anInputId);
335               }
336             }
337           }
338
339           anOutputCellData->AddArray(aDataArray);
340           aDataArray->Delete();
341         }
342         return true;
343       }
344     }
345
346     return false;
347   }
348
349
350   //---------------------------------------------------------------
351 }
352
353
354 namespace VISU
355 {
356   //---------------------------------------------------------------
357   TAppendFilterHelper
358   ::TAppendFilterHelper(vtkObject* theParent):
359     myIsMergingInputs(false),
360     myIsMappingInputs(false),
361     myParent(*theParent)
362   {}
363
364
365   //---------------------------------------------------------------
366   void
367   TAppendFilterHelper
368   ::SetSharedPointSet(vtkPointSet* thePointSet)
369   {
370     if(GetSharedPointSet() == thePointSet)
371       return;
372     
373     mySharedPointSet = thePointSet;
374     
375     myParent.Modified();
376   }
377
378
379   //---------------------------------------------------------------
380   vtkPointSet*
381   TAppendFilterHelper
382   ::GetSharedPointSet()
383   {
384     return mySharedPointSet.GetPointer();
385   }
386   
387
388   //---------------------------------------------------------------
389   void
390   TAppendFilterHelper
391   ::SetMappingInputs(bool theIsMappingInputs)
392   {
393     if(myIsMappingInputs == theIsMappingInputs)
394       return;
395     
396     myIsMappingInputs = theIsMappingInputs;
397     myParent.Modified();
398   }
399   
400   
401   //---------------------------------------------------------------
402   bool
403   TAppendFilterHelper
404   ::IsMappingInputs()
405   {
406     return myIsMappingInputs;
407   }
408   
409
410   //---------------------------------------------------------------
411   void
412   TAppendFilterHelper
413   ::SetMergingInputs(bool theIsMergingInputs)
414   {
415     if(myIsMergingInputs == theIsMergingInputs)
416       return;
417     
418     myIsMergingInputs = theIsMergingInputs;
419     myParent.Modified();
420   }
421   
422   
423   //---------------------------------------------------------------
424   bool
425   TAppendFilterHelper
426   ::IsMergingInputs()
427   {
428     return myIsMergingInputs;
429   }
430   
431
432   //---------------------------------------------------------------
433   bool
434   UnstructuredGridRequestData(vtkInformationVector **theInputVector,
435                               vtkIdType theNumberOfInputConnections,
436                               vtkInformationVector *theOutputVector,
437                               vtkPointSet* theSharedPointSet,
438                               bool theIsMergingInputs,
439                               bool theIsMappingInputs)
440   {
441     return RequestData<vtkUnstructuredGrid>(theInputVector,
442                                             theNumberOfInputConnections,
443                                             theOutputVector,
444                                             theSharedPointSet,
445                                             theIsMergingInputs,
446                                             theIsMappingInputs);
447   }
448
449
450   //---------------------------------------------------------------
451   bool
452   PolyDataRequestData(vtkInformationVector **theInputVector,
453                       vtkIdType theNumberOfInputConnections,
454                       vtkInformationVector *theOutputVector,
455                       vtkPointSet* theSharedPointSet,
456                       bool theIsMergingInputs,
457                       bool theIsMappingInputs)
458   {
459     return RequestData<vtkPolyData>(theInputVector,
460                                     theNumberOfInputConnections,
461                                     theOutputVector,
462                                     theSharedPointSet,
463                                     theIsMergingInputs,
464                                     theIsMappingInputs);
465   }
466
467
468   //---------------------------------------------------------------
469 }