]> SALOME platform Git repositories - modules/visu.git/blob - src/CONVERTOR/VISU_MergeFilterUtilities.cxx
Salome HOME
adcf7fb0aedeae472234ae466bdde3dd7216b7e3
[modules/visu.git] / src / CONVERTOR / VISU_MergeFilterUtilities.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 VTKViewer : build VTK viewer into Salome desktop
23 //  File   : 
24 //  Author : 
25 //  Module : SALOME
26 //  $Header$
27 //
28 #include "VISU_MergeFilterUtilities.hxx"
29
30 #include <vtkCellData.h>
31 #include <vtkObjectFactory.h>
32 #include <vtkPointData.h>
33 #include <vtkPolyData.h>
34 #include <vtkRectilinearGrid.h>
35 #include <vtkStructuredGrid.h>
36 #include <vtkStructuredPoints.h>
37 #include <vtkUnstructuredGrid.h>
38
39 #include <vtkIdList.h>
40 #include <vtkCell.h>
41
42 #include <algorithm>
43 #include <vector>
44 #include <set>
45 #include <map>
46
47 #if !defined(VTK_XVERSION)
48 #define VTK_XVERSION (VTK_MAJOR_VERSION<<16)+(VTK_MINOR_VERSION<<8)+(VTK_BUILD_VERSION)
49 #endif
50
51 namespace
52 {
53
54   using namespace VISU;
55   
56   //---------------------------------------------------------------
57   template<class TDataSet>
58   
59   void 
60   CopyDataOnCells(TDataSet *theInput,
61                   vtkIntArray *theGeometryCellMapper,
62                   vtkIntArray *theDataCellMapper,
63                   vtkDataSet* theScalarsDataSet,
64                   vtkDataSet* theVectorsDataSet,
65                   vtkDataSet* theNormalsDataSet,
66                   vtkDataSet* theTCoordsDataSet,
67                   vtkDataSet* theTensorsDataSet,
68                   VISU::TFieldList* theFieldList,
69                   TDataSet *theOutput)
70   {
71     if(IsDifferent(theGeometryCellMapper, theDataCellMapper)){
72       TObjectIdArray anIntersection;
73       GetIntersection(theGeometryCellMapper,
74                       theDataCellMapper,
75                       anIntersection);
76     
77       TObjectId2TupleIdMap aGeomObjectId2TupleIdMap;
78       GetObjectId2TupleIdMap(theGeometryCellMapper, aGeomObjectId2TupleIdMap);
79       
80       TObjectId2TupleIdMap aDataObjectId2TupleIdMap;
81       GetObjectId2TupleIdMap(theDataCellMapper, aDataObjectId2TupleIdMap);
82       
83       vtkCellData *aCellData = theScalarsDataSet->GetCellData();
84       vtkCellData *anOutputCellData = theOutput->GetCellData();
85       anOutputCellData->CopyAllocate(aCellData);
86
87       if(theVectorsDataSet && theVectorsDataSet != theScalarsDataSet)
88         anOutputCellData->CopyVectorsOff();
89
90       vtkIdType aNbTuples = anIntersection.size();
91       theOutput->Allocate(aNbTuples);
92       vtkIdList *aCellIds = vtkIdList::New();
93       for(int aTupleId = 0; aTupleId < aNbTuples; aTupleId++){
94         TObjectId& anObjectId = anIntersection[aTupleId];
95         vtkIdType aCellId = aGeomObjectId2TupleIdMap[anObjectId];
96         vtkCell *aCell = theInput->GetCell(aCellId);
97         aCellIds->Reset();
98         vtkIdType aNbPointIds = aCell->PointIds->GetNumberOfIds();
99         for(vtkIdType anId = 0; anId < aNbPointIds; anId++)
100           aCellIds->InsertNextId(aCell->GetPointIds()->GetId(anId));
101         vtkIdType aCellType = theInput->GetCellType(aCellId);
102         vtkIdType aNewCellId = theOutput->InsertNextCell(aCellType, aCellIds);
103         vtkIdType aDataCellId = aDataObjectId2TupleIdMap[anObjectId];
104         anOutputCellData->CopyData(aCellData, aDataCellId, aNewCellId);
105       }
106       aCellIds->Delete();
107
108       theOutput->SetPoints(theInput->GetPoints());
109       
110     }else{
111       theOutput->CopyStructure(theInput);
112       theOutput->GetCellData()->ShallowCopy(theScalarsDataSet->GetCellData());
113     }
114     theOutput->GetPointData()->ShallowCopy(theInput->GetPointData());
115     
116     //If need, copy vectors data.
117     if(theVectorsDataSet && theVectorsDataSet != theScalarsDataSet){
118       bool isVectorsOnCells = theVectorsDataSet->GetCellData()->GetVectors() != NULL;
119       bool isVectorsDataOnPoints = theVectorsDataSet->GetPointData()->GetVectors() != NULL;
120       if(isVectorsOnCells) {
121         CopyVectorsOnCells(theVectorsDataSet,theOutput);
122       }
123       else if(isVectorsDataOnPoints){
124         CopyVectorsOnPoints(theVectorsDataSet,theOutput);
125       }
126     }
127   }
128   
129   void CopyVectorsOnCells(vtkDataSet *theVectorsDataSet,
130                           vtkDataSet *theOutput)
131   {
132     vtkDataArray *anInputVectors = theVectorsDataSet->GetCellData()->GetVectors();
133     vtkDataArray *anOutputVectors = vtkDataArray::CreateDataArray(anInputVectors->GetDataType());
134     
135     //Clear output vector data
136     theOutput->GetCellData()->SetVectors(NULL);
137     
138     //Copy vectors data
139     vtkIntArray* anOutputIDMapper = GetIDMapper(theOutput,
140                                                 TGetCellData(),
141                                                 "VISU_CELLS_MAPPER");
142     
143     vtkIntArray* anInputIDMapper = GetIDMapper(theVectorsDataSet,
144                                                TGetCellData(),
145                                                "VISU_CELLS_MAPPER");
146     
147     TObjectIdArray anIntersection;
148     GetIntersection(anOutputIDMapper,
149                     anInputIDMapper,
150                     anIntersection);
151
152     vtkIdType aNbTuples = anIntersection.size();
153     anOutputVectors->SetNumberOfComponents(anInputVectors->GetNumberOfComponents());
154     anOutputVectors->SetNumberOfTuples(aNbTuples);
155     theOutput->GetCellData()->SetVectors(anOutputVectors);
156     anOutputVectors->Delete();
157     
158     TObjectId2TupleIdMap anOutputObjectId2TupleIdMap;
159     GetObjectId2TupleIdMap(anOutputIDMapper, anOutputObjectId2TupleIdMap);
160     
161     TObjectId2TupleIdMap anInputObjectId2TupleIdMap;
162     GetObjectId2TupleIdMap(anInputIDMapper, anInputObjectId2TupleIdMap);
163
164     for(vtkIdType iTupleId = 0; iTupleId < aNbTuples; iTupleId++ ){
165       TObjectId &anObjectId = anIntersection[iTupleId];
166       vtkIdType anOutputCellId  = anOutputObjectId2TupleIdMap[anObjectId];
167       vtkIdType anInputCellId = anInputObjectId2TupleIdMap[anObjectId];
168       anOutputVectors->SetTuple(anOutputCellId,anInputVectors->GetTuple(anInputCellId));
169     }
170   }
171   
172   //---------------------------------------------------------------
173   template<class TDataSet>
174   void 
175   CopyDataOnPoints(TDataSet *theInput,
176                    vtkIntArray *theGeometryPointMapper,
177                    vtkIntArray *theDataPointMapper,
178                    vtkDataSet* theScalarsDataSet,
179                    vtkDataSet* theVectorsDataSet,
180                    vtkDataSet* theNormalsDataSet,
181                    vtkDataSet* theTCoordsDataSet,
182                    vtkDataSet* theTensorsDataSet,
183                    VISU::TFieldList* theFieldList,
184                    TDataSet *theOutput)
185   {
186     if(IsDifferent(theGeometryPointMapper, theDataPointMapper)){
187       TObjectId2TupleIdMap aDataObjectId2PointIdMap;
188       GetObjectId2TupleIdMap(theDataPointMapper, aDataObjectId2PointIdMap);
189
190       vtkCellData *aCellData = theInput->GetCellData();
191       vtkCellData *anOutputCellData = theOutput->GetCellData();
192       anOutputCellData->CopyAllocate(aCellData);
193
194       vtkIdList *aCellIds = vtkIdList::New();
195       int aNbCells = theInput->GetNumberOfCells();
196       theOutput->Allocate(aNbCells);
197       for(int aCellId = 0; aCellId < aNbCells; aCellId++){
198         aCellIds->Reset();
199         vtkCell *aCell = theInput->GetCell(aCellId);
200         vtkIdType aNbPointIds = aCell->PointIds->GetNumberOfIds();
201         for(vtkIdType anId = 0; anId < aNbPointIds; anId++){
202           vtkIdType aPointId = aCell->GetPointIds()->GetId(anId);
203           int* aPointer = theGeometryPointMapper->GetPointer(aPointId * 2);
204           TCellId aCellId = *aPointer;
205           TEntityId anEntityId = *(aPointer + 1);
206           TObjectId anObjectId(aCellId, anEntityId);
207           TObjectId2TupleIdMap::iterator anIter = aDataObjectId2PointIdMap.find(anObjectId);
208           if(anIter != aDataObjectId2PointIdMap.end()){
209             aPointId = anIter->second;
210             aCellIds->InsertNextId(aPointId);
211           }else
212             goto PASS_INSERT_NEXT_CELL;
213         }
214         {
215           vtkIdType aCellType = theInput->GetCellType(aCellId);
216           vtkIdType aNewCellId = theOutput->InsertNextCell(aCellType, aCellIds);
217           anOutputCellData->CopyData(aCellData, aCellId, aNewCellId);
218         }
219       PASS_INSERT_NEXT_CELL:
220         continue;
221       }
222       aCellIds->Delete();
223       
224       // Copy geometry points
225       // 1. Create vtkPoints instance of the same data type
226       vtkPointSet* aScalarsDataSet = dynamic_cast<vtkPointSet*>(theScalarsDataSet);
227       vtkPoints* anGeometryPoints = theInput->GetPoints();
228       vtkPoints* aDataPoints = aScalarsDataSet->GetPoints();
229       vtkPoints* anOutputPoints = vtkPoints::New(aDataPoints->GetDataType());
230       theOutput->SetPoints(anOutputPoints);
231       anOutputPoints->Delete();
232
233       // 2. Perform mapping of geometry points
234       TObjectId2TupleIdMap aGeomObjectId2TupleIdMap;
235       GetObjectId2TupleIdMap(theGeometryPointMapper, aGeomObjectId2TupleIdMap);
236
237       // 3. Loop over all data points
238       int aNbDataPoints = theDataPointMapper->GetNumberOfTuples();
239       anOutputPoints->SetNumberOfPoints(aNbDataPoints);
240       for(int aPointId = 0; aPointId < aNbDataPoints; aPointId++){
241         int* aPointer = theDataPointMapper->GetPointer(aPointId * 2);
242         TCellId aCellId = *aPointer;
243         TEntityId anEntityId = *(aPointer + 1);
244         TObjectId anObjectId(aCellId, anEntityId);
245         TObjectId2TupleIdMap::iterator anIter = aGeomObjectId2TupleIdMap.find(anObjectId);
246         if(anIter != aGeomObjectId2TupleIdMap.end()){
247           // If the point exists in the geometry put it to output
248           int aGeometryPointId = anIter->second;
249           vtkFloatingPointType aCoords[3];
250           anGeometryPoints->GetPoint(aGeometryPointId, aCoords);
251           anOutputPoints->SetPoint(aPointId, aCoords);
252         }else{
253           // If no, the point from data should be used
254           vtkFloatingPointType aCoords[3];
255           aDataPoints->GetPoint(aPointId, aCoords);
256           anOutputPoints->SetPoint(aPointId, aCoords);
257         }
258       }
259     }else{
260       theOutput->CopyStructure(theInput);
261       theOutput->GetCellData()->ShallowCopy(theInput->GetCellData());
262     }
263     theOutput->GetPointData()->ShallowCopy(theScalarsDataSet->GetPointData());
264     
265     //If need, copy vectors data.
266     if(theVectorsDataSet && theVectorsDataSet != theScalarsDataSet){
267       bool isVectorsOnCells = theVectorsDataSet->GetCellData()->GetVectors() != NULL;
268       bool isVectorsDataOnPoints = theVectorsDataSet->GetPointData()->GetVectors() != NULL;
269
270       //Merge cells if need
271       //rnv
272       if(!IsDifferent(theGeometryPointMapper, theDataPointMapper)){
273         vtkIntArray* theGeometryCellMapper = GetIDMapper(theVectorsDataSet,
274                                                          TGetCellData(),
275                                                          "VISU_CELLS_MAPPER");
276         
277         vtkIntArray* theDataCellMapper = GetIDMapper(theScalarsDataSet,
278                                                      TGetCellData(),
279                                                      "VISU_CELLS_MAPPER");
280         
281           
282         if(IsDifferent(theGeometryCellMapper, theDataCellMapper)){
283           TObjectIdArray anIntersection;
284           
285           GetIntersection(theGeometryCellMapper,
286                           theDataCellMapper,
287                           anIntersection);
288     
289           TObjectId2TupleIdMap aGeomObjectId2TupleIdMap;
290           GetObjectId2TupleIdMap(theGeometryCellMapper, aGeomObjectId2TupleIdMap);
291       
292           TObjectId2TupleIdMap aDataObjectId2TupleIdMap;
293           GetObjectId2TupleIdMap(theDataCellMapper, aDataObjectId2TupleIdMap);
294       
295           vtkCellData *aCellData = theScalarsDataSet->GetCellData();
296           vtkCellData *anOutputCellData = theOutput->GetCellData();
297           anOutputCellData->CopyAllocate(aCellData);
298
299           vtkIdType aNbTuples = anIntersection.size();
300           theOutput->Allocate(aNbTuples);
301           vtkIdList *aCellIds = vtkIdList::New();
302           for(int aTupleId = 0; aTupleId < aNbTuples; aTupleId++){
303             TObjectId& anObjectId = anIntersection[aTupleId];
304             vtkIdType aCellId = aGeomObjectId2TupleIdMap[anObjectId];
305             vtkCell *aCell = theInput->GetCell(aCellId);
306             aCellIds->Reset();
307             vtkIdType aNbPointIds = aCell->PointIds->GetNumberOfIds();
308             for(vtkIdType anId = 0; anId < aNbPointIds; anId++)
309               aCellIds->InsertNextId(aCell->GetPointIds()->GetId(anId));
310             vtkIdType aCellType = theInput->GetCellType(aCellId);
311             vtkIdType aNewCellId = theOutput->InsertNextCell(aCellType, aCellIds);
312             vtkIdType aDataCellId = aDataObjectId2TupleIdMap[anObjectId];
313             anOutputCellData->CopyData(aCellData, aDataCellId, aNewCellId);
314           }
315           aCellIds->Delete();
316           
317         }
318       }
319         
320       if(isVectorsOnCells) {
321         CopyVectorsOnCells(theVectorsDataSet,theOutput);
322       }
323       else if(isVectorsDataOnPoints){
324         CopyVectorsOnPoints(theVectorsDataSet,theOutput);
325       }
326     }
327   }
328   
329   void CopyVectorsOnPoints(vtkDataSet *theVectorsDataSet,
330                           vtkDataSet *theOutput)
331   {
332     vtkDataArray *anInputVectors = theVectorsDataSet->GetPointData()->GetVectors();
333     
334     //Clear output vector data
335     theOutput->GetPointData()->SetVectors(NULL);
336     
337     vtkDataArray *anOutputVectors = vtkDataArray::CreateDataArray(anInputVectors->GetDataType());
338     
339     //Copy vectors data
340     vtkIntArray* anOutputIDMapper = GetIDMapper(theOutput,
341                                                 TGetPointData(),
342                                                 "VISU_POINTS_MAPPER");
343     
344     vtkIntArray* anInputIDMapper = GetIDMapper(theVectorsDataSet,
345                                                TGetPointData(),
346                                                "VISU_POINTS_MAPPER");
347     TObjectIdArray anIntersection;
348
349     GetIntersection(anOutputIDMapper,
350                     anInputIDMapper,
351                     anIntersection);
352     
353     vtkIdType aNbTuples = anIntersection.size();  
354     anOutputVectors->SetNumberOfComponents(anInputVectors->GetNumberOfComponents());
355     anOutputVectors->SetNumberOfTuples(aNbTuples);
356     
357     
358
359     TObjectId2TupleIdMap anOutputObjectId2TupleIdMap;
360     GetObjectId2TupleIdMap(anOutputIDMapper, anOutputObjectId2TupleIdMap);
361     
362     TObjectId2TupleIdMap anInputObjectId2TupleIdMap;
363     GetObjectId2TupleIdMap(anInputIDMapper, anInputObjectId2TupleIdMap);
364     
365     for(vtkIdType iTupleId = 0; iTupleId < aNbTuples; iTupleId++ ){
366       TObjectId& anObjectId = anIntersection[iTupleId];
367       vtkIdType anOutputPointId  = anOutputObjectId2TupleIdMap[anObjectId];
368       vtkIdType anInputPointId = anInputObjectId2TupleIdMap[anObjectId];
369       anOutputVectors->SetTuple(anOutputPointId,anInputVectors->GetTuple(anInputPointId));
370     }
371     
372     theOutput->GetPointData()->SetVectors(anOutputVectors);
373     anOutputVectors->Delete();
374   }
375   
376
377   //---------------------------------------------------------------
378   typedef vtkDataArray* (vtkDataSetAttributes::* TGetAttribute)();
379 #if (VTK_XVERSION < 0x050100)
380   typedef int (vtkDataSetAttributes::* TSetAttribute)(vtkDataArray*);
381 #else
382   typedef int (vtkDataSetAttributes::* TSetAttribute)(vtkAbstractArray*);
383   typedef int (vtkDataSetAttributes::* TSetDataAttribute)(vtkDataArray*);
384 #endif
385
386
387 #if (VTK_XVERSION >= 0x050100)
388   inline
389   void
390   CopyArray(vtkDataArray* theDataArray,
391             vtkDataSetAttributes* theOutput, 
392             TSetDataAttribute theSetAttribute,
393             vtkIdType theFixedNbTuples)
394   {
395     if(theDataArray){
396       vtkIdType aNbTuples = theDataArray->GetNumberOfTuples();
397       if(theFixedNbTuples == aNbTuples)
398         (theOutput->*theSetAttribute)(theDataArray);
399     }
400   }
401 #endif
402
403   inline
404   void
405   CopyArray(vtkDataArray* theDataArray,
406             vtkDataSetAttributes* theOutput, 
407             TSetAttribute theSetAttribute,
408             vtkIdType theFixedNbTuples)
409   {
410     if(theDataArray){
411       vtkIdType aNbTuples = theDataArray->GetNumberOfTuples();
412       if(theFixedNbTuples == aNbTuples)
413         (theOutput->*theSetAttribute)(theDataArray);
414     }
415   }
416
417
418   //---------------------------------------------------------------
419   inline
420   void
421   CopyAttribute(vtkDataSetAttributes* theInput, 
422                 TGetAttribute theGetAttribute,
423                 vtkDataSetAttributes* theOutput, 
424 #if (VTK_XVERSION < 0x050100)
425                 TSetAttribute theSetAttribute,
426 #else
427                 TSetDataAttribute theSetAttribute,
428 #endif
429                 vtkIdType theFixedNbTuples)
430   {
431     CopyArray((theInput->*theGetAttribute)(),
432               theOutput, theSetAttribute,
433               theFixedNbTuples);
434   }
435
436
437   //---------------------------------------------------------------
438   inline
439   void
440   CopyDataSetAttribute(vtkDataSet *theInput,
441                        TGetAttribute theGetAttribute,
442                        vtkDataSet *theOutput, 
443 #if (VTK_XVERSION < 0x050100)
444                        TSetAttribute theSetAttribute,
445 #else
446                        TSetDataAttribute theSetAttribute,
447 #endif
448                        vtkIdType theNbPoints, 
449                        vtkIdType theNbCells)
450   {
451     CopyAttribute(theInput->GetPointData(), 
452                   theGetAttribute,
453                   theOutput->GetPointData(), 
454                   theSetAttribute,
455                   theNbPoints);
456     CopyAttribute(theInput->GetCellData(), 
457                   theGetAttribute,
458                   theOutput->GetCellData(), 
459                   theSetAttribute,
460                   theNbCells);
461   }
462
463
464   //---------------------------------------------------------------
465   inline
466   void
467   CopyField(vtkDataSetAttributes* theInput, 
468             const char* theFieldName, 
469             vtkDataSetAttributes* theOutput,
470             vtkIdType theFixedNbTuples)
471   {
472     CopyArray(theInput->GetArray(theFieldName),
473               theOutput,
474 #if (VTK_XVERSION < 0x050100)
475         &vtkDataSetAttributes::AddArray,
476 #else
477               &vtkFieldData::AddArray,
478 #endif
479               theFixedNbTuples);
480   }
481
482
483   //---------------------------------------------------------------
484   inline
485   void
486   CopyDataSetField(vtkDataSet* theInput, 
487                    const char* theFieldName, 
488                    vtkDataSet* theOutput,
489                    vtkIdType theNbPoints, 
490                    vtkIdType theNbCells)
491   {
492     if(theInput){
493       CopyField(theInput->GetPointData(), 
494                 theFieldName, 
495                 theOutput->GetPointData(), 
496                 theNbPoints);
497       CopyField(theInput->GetCellData(), 
498                 theFieldName, 
499                 theOutput->GetCellData(), 
500                 theNbCells);
501     }
502   }
503
504   //---------------------------------------------------------------
505   void
506   BasicExecute(vtkDataSet *theInput,
507                vtkDataSet* theScalarsDataSet,
508                vtkDataSet* theVectorsDataSet,
509                vtkDataSet* theNormalsDataSet,
510                vtkDataSet* theTCoordsDataSet,
511                vtkDataSet* theTensorsDataSet,
512                VISU::TFieldList* theFieldList,
513                vtkDataSet *theOutput)
514   {
515     theOutput->CopyStructure(theInput);
516
517     vtkIdType aNbPoints = theInput->GetNumberOfPoints();
518     vtkIdType aNbCells = theInput->GetNumberOfCells();
519   
520     // merge data only if it is consistent
521     if(theScalarsDataSet)
522       CopyDataSetAttribute(theScalarsDataSet, 
523                            &vtkDataSetAttributes::GetScalars,
524                            theOutput,
525                            &vtkDataSetAttributes::SetScalars,
526                            aNbPoints,
527                            aNbCells);
528
529     if(theVectorsDataSet)
530       CopyDataSetAttribute(theVectorsDataSet, 
531                            &vtkDataSetAttributes::GetVectors,
532                            theOutput, 
533                            &vtkDataSetAttributes::SetVectors,
534                            aNbPoints, 
535                            aNbCells);
536
537     if(theNormalsDataSet)
538       CopyDataSetAttribute(theNormalsDataSet,
539                            &vtkDataSetAttributes::GetNormals,
540                            theOutput,
541                            &vtkDataSetAttributes::SetNormals,
542                            aNbPoints,
543                            aNbCells);
544
545     if(theTCoordsDataSet)
546       CopyDataSetAttribute(theTCoordsDataSet,
547                            &vtkDataSetAttributes::GetTCoords,
548                            theOutput,
549                            &vtkDataSetAttributes::SetTCoords,
550                            aNbPoints, 
551                            aNbCells);
552     
553     if(theTensorsDataSet)
554       CopyDataSetAttribute(theTensorsDataSet, 
555                            &vtkDataSetAttributes::GetTensors,
556                            theOutput, 
557                            &vtkDataSetAttributes::SetTensors,
558                            aNbPoints, 
559                            aNbCells);
560
561     VISU::TFieldListIterator anIter(theFieldList);
562     for(anIter.Begin(); !anIter.End() ; anIter.Next()){
563       vtkDataSet *aDataSet = anIter.Get()->Ptr;
564       const char* aFieldName = anIter.Get()->GetName();
565       CopyDataSetField(aDataSet, 
566                        aFieldName, 
567                        theOutput, 
568                        aNbPoints, 
569                        aNbCells);
570     }
571   }
572
573
574   //---------------------------------------------------------------
575   template<class TDataSet>
576   bool 
577   Execute(TDataSet *theInput,
578           vtkDataSet* theScalarsDataSet,
579           vtkDataSet* theVectorsDataSet,
580           vtkDataSet* theNormalsDataSet,
581           vtkDataSet* theTCoordsDataSet,
582           vtkDataSet* theTensorsDataSet,
583           VISU::TFieldList* theFieldList,
584           bool theIsMergingInputs,
585           TDataSet *theOutput)
586   {
587     if(theIsMergingInputs){
588       vtkCellData *aCellData = theInput->GetCellData();
589       if(vtkDataArray *aCellMapper = aCellData->GetArray("VISU_CELLS_MAPPER")){
590         bool anIsDataOnCells = false;
591         if(vtkDataSet* aDataSet = theScalarsDataSet)
592           if(vtkCellData* aCellData = aDataSet->GetCellData())
593             anIsDataOnCells = (aCellData->GetArray("VISU_FIELD") != NULL);
594         if(anIsDataOnCells){
595           vtkIntArray *aGeometryCellMapper = dynamic_cast<vtkIntArray*>(aCellMapper);
596           vtkIntArray* aDataCellMapper = GetIDMapper(theFieldList,
597                                                      TGetCellData(),
598                                                      "VISU_CELLS_MAPPER");
599           CopyDataOnCells(theInput,
600                           aGeometryCellMapper,
601                           aDataCellMapper,
602                           theScalarsDataSet,
603                           theVectorsDataSet,
604                           theNormalsDataSet,
605                           theTCoordsDataSet,
606                           theTensorsDataSet,
607                           theFieldList,
608                           theOutput);
609         }else{
610           vtkPointData *aPointData = theInput->GetPointData();
611           vtkDataArray *aPointMapper = aPointData->GetArray("VISU_POINTS_MAPPER");
612           vtkIntArray *aGeometryPointMapper = dynamic_cast<vtkIntArray*>(aPointMapper);
613           vtkIntArray* aDataPointMapper = GetIDMapper(theFieldList,
614                                                       TGetPointData(),
615                                                       "VISU_POINTS_MAPPER");
616           CopyDataOnPoints(theInput,
617                            aGeometryPointMapper,
618                            aDataPointMapper,
619                            theScalarsDataSet,
620                            theVectorsDataSet,
621                            theNormalsDataSet,
622                            theTCoordsDataSet,
623                            theTensorsDataSet,
624                            theFieldList,
625                            theOutput);
626         }
627       }
628     }else{
629       BasicExecute(theInput,
630                    theScalarsDataSet,
631                    theVectorsDataSet,
632                    theNormalsDataSet,
633                    theTCoordsDataSet,
634                    theTensorsDataSet,
635                    theFieldList,
636                    theOutput);
637     }
638     return true;
639   }
640 }
641
642 namespace VISU
643 {
644
645   //---------------------------------------------------------------
646   void
647   GetObjectIdSet(vtkIntArray *theArray, 
648                  TObjectIdSet& theObjectIdSet)
649   {
650     theObjectIdSet.clear();
651     int aMaxId = theArray->GetMaxId();
652     int* aPointer = theArray->GetPointer(0);
653     int* anEndPointer = theArray->GetPointer(aMaxId + 1);
654     for(; aPointer != anEndPointer; aPointer += 2){
655       TCellId aCellId = *aPointer;
656       TEntityId anEntityId = *(aPointer + 1);
657       TObjectId anObjectId(aCellId, anEntityId);
658       theObjectIdSet.insert(anObjectId);
659     }
660   }
661
662
663   //---------------------------------------------------------------
664   void
665   GetObjectId2TupleIdMap(vtkIntArray *theArray, 
666                         TObjectId2TupleIdMap& theObjectId2TupleIdMap)
667   {
668     theObjectId2TupleIdMap.clear();
669     int* aPointer = theArray->GetPointer(0);
670     int aNbTuples = theArray->GetNumberOfTuples();
671     for(vtkIdType aTupleId = 0; aTupleId < aNbTuples; aTupleId++, aPointer += 2){
672       TCellId aCellId = *aPointer;
673       TEntityId anEntityId = *(aPointer + 1);
674       TObjectId anObjectId(aCellId, anEntityId);
675       theObjectId2TupleIdMap[anObjectId] = aTupleId;
676     }
677   }
678
679   //---------------------------------------------------------------
680   void 
681   GetObjectId2TupleGaussIdArray(vtkIntArray *theArray, 
682                                 TObjectId2TupleGaussIdMap& theObjectId2TupleGaussIdMap)
683   {
684     theObjectId2TupleGaussIdMap.clear();
685     int* aPointer = theArray->GetPointer(0);
686     int aNbTuples = theArray->GetNumberOfTuples();
687     for(vtkIdType aTupleId = 0; aTupleId < aNbTuples; aTupleId++){
688       int aCellId = *aPointer;
689       TObjectId2TupleGaussIdMap::iterator it = theObjectId2TupleGaussIdMap.find(aCellId);
690       if(it == theObjectId2TupleGaussIdMap.end()){
691         TCellIdArray anIdArray;
692         anIdArray.push_back(aTupleId);
693         theObjectId2TupleGaussIdMap.insert(make_pair(aCellId,anIdArray));
694       }
695       else{
696         (*it).second.push_back(aTupleId);
697       }
698       aPointer += 2;
699     }
700   }
701   
702   //---------------------------------------------------------------
703   template<class TGetFieldData>
704   vtkIntArray*
705   GetIDMapper(VISU::TFieldList* theFieldList,
706               TGetFieldData theGetFieldData,
707               const char* theFieldName)
708   {
709     VISU::TFieldListIterator anIter(theFieldList);
710     for(anIter.Begin(); !anIter.End() ; anIter.Next()){
711       const char* aFieldName = anIter.Get()->GetName();
712       if(strcmp(aFieldName, theFieldName) == 0){
713         vtkDataSet* aDataSet = anIter.Get()->Ptr;
714         vtkFieldData *aFieldData = theGetFieldData(aDataSet);
715         vtkDataArray *anIDMapper = aFieldData->GetArray(theFieldName);
716         return dynamic_cast<vtkIntArray*>(anIDMapper);
717       }
718     }
719     return NULL;
720   }
721
722
723   //---------------------------------------------------------------
724   template<class TGetFieldData>
725   vtkIntArray*
726   GetIDMapper(vtkDataSet* theIDMapperDataSet,
727               TGetFieldData theGetFieldData,
728               const char* theFieldName)
729   {
730     vtkFieldData *aFieldData = theGetFieldData(theIDMapperDataSet);
731     vtkDataArray *anIDMapper = aFieldData->GetArray(theFieldName);
732     return dynamic_cast<vtkIntArray*>(anIDMapper);
733   }
734
735
736   //---------------------------------------------------------------
737   bool
738   IsDifferent(vtkIntArray *theFirstIDMapper,
739               vtkIntArray *theSecondIDMapper)
740   {
741     vtkIdType aFirstNbTuples = theFirstIDMapper->GetNumberOfTuples();
742     vtkIdType aSecondNbTuples = theSecondIDMapper->GetNumberOfTuples();
743     if(aFirstNbTuples != aSecondNbTuples)
744       return true;
745
746     int aMaxId = theFirstIDMapper->GetMaxId();
747     int* aFirstPointer = theFirstIDMapper->GetPointer(0);
748     int* aSecondPointer = theSecondIDMapper->GetPointer(0);
749     for(int anId = 0; anId <= aMaxId; anId++){
750       if(*aFirstPointer++ != *aSecondPointer++)
751         return true;
752     }
753     
754     return false;
755   }
756
757
758   //---------------------------------------------------------------
759   void
760   GetIntersection(vtkIntArray *theFirstIDMapper,
761                   vtkIntArray *theSecondIDMapper,
762                   TObjectIdArray& theResult)
763   {
764     TObjectIdSet aFirstObjectIdSet;
765     GetObjectIdSet(theFirstIDMapper, aFirstObjectIdSet);
766     
767     TObjectIdSet aSecondObjectIdSet;
768     GetObjectIdSet(theSecondIDMapper, aSecondObjectIdSet);
769
770     size_t aMaxLength = std::max(aFirstObjectIdSet.size(), aSecondObjectIdSet.size());
771     theResult.resize(aMaxLength);
772     TObjectIdArray::iterator anArrayIter = theResult.begin();
773     anArrayIter = std::set_intersection(aFirstObjectIdSet.begin(),
774                                         aFirstObjectIdSet.end(),
775                                         aSecondObjectIdSet.begin(),
776                                         aSecondObjectIdSet.end(),
777                                         anArrayIter);
778     theResult.erase(anArrayIter, theResult.end());
779   }
780
781   //---------------------------------------------------------------
782   bool
783   Execute(vtkUnstructuredGrid *theInput,
784           vtkUnstructuredGrid *theOutput,
785           vtkDataSet* theScalarsDataSet,
786           vtkDataSet* theVectorsDataSet,
787           vtkDataSet* theNormalsDataSet,
788           vtkDataSet* theTCoordsDataSet,
789           vtkDataSet* theTensorsDataSet,
790           TFieldList* theFieldList,
791           bool theIsMergingInputs)
792   {
793     return ::Execute(theInput,
794                      theScalarsDataSet,
795                      theVectorsDataSet,
796                      theNormalsDataSet,
797                      theTCoordsDataSet,
798                      theTensorsDataSet,
799                      theFieldList,
800                      theIsMergingInputs,
801                      theOutput);
802   }
803
804
805   //---------------------------------------------------------------
806   bool
807   Execute(vtkPolyData *theInput,
808           vtkPolyData *theOutput,
809           vtkDataSet* theScalarsDataSet,
810           vtkDataSet* theVectorsDataSet,
811           vtkDataSet* theNormalsDataSet,
812           vtkDataSet* theTCoordsDataSet,
813           vtkDataSet* theTensorsDataSet,
814           TFieldList* theFieldList,
815           bool theIsMergingInputs)
816   {
817     return ::Execute(theInput,
818                      theScalarsDataSet,
819                      theVectorsDataSet,
820                      theNormalsDataSet,
821                      theTCoordsDataSet,
822                      theTensorsDataSet,
823                      theFieldList,
824                      theIsMergingInputs,
825                      theOutput);
826   }
827 }