]> SALOME platform Git repositories - modules/visu.git/blob - src/PIPELINE/VISU_ElnoDisassembleFilter.cxx
Salome HOME
38957e4bb786160d648a48b3765b822c9b3f44b1
[modules/visu.git] / src / PIPELINE / VISU_ElnoDisassembleFilter.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 #include "VISU_ElnoDisassembleFilter.hxx"
23 #include "VISU_PipeLineUtils.hxx"
24 #include "VISU_ElnoMeshValue.hxx"
25
26 #include <vtkCellData.h>
27 #include <vtkInformation.h>
28 #include <vtkInformationVector.h>
29 #include <vtkObjectFactory.h>
30 #include <vtkPointData.h>
31 #include <vtkUnstructuredGrid.h>
32 #include <vtkPoints.h>
33 #include <vtkCellArray.h>
34
35
36 //----------------------------------------------------------------------------
37 vtkStandardNewMacro( VISU_ElnoDisassembleFilter );
38
39
40 //----------------------------------------------------------------------------
41 VISU_ElnoDisassembleFilter::VISU_ElnoDisassembleFilter()
42 {
43   this->SetInputArrayToProcess( 0, // idx
44                                 0, // port
45                                 0, // connection
46                                 vtkDataObject::FIELD_ASSOCIATION_CELLS, // field association
47                                 "ELNO_FIELD" ); // name
48
49   this->SetInputArrayToProcess( 1, // idx
50                                 0, // port
51                                 0, // connection
52                                 vtkDataObject::FIELD_ASSOCIATION_CELLS, // field association
53                                 "ELNO_COMPONENT_MAPPER" ); // name
54
55   this->myShrinkFactor = -0.999;
56 }
57
58
59 //----------------------------------------------------------------------------
60 VISU_ElnoDisassembleFilter::~VISU_ElnoDisassembleFilter()
61 {}
62
63
64 //----------------------------------------------------------------------------
65 void VISU_ElnoDisassembleFilter::SetShrinkFactor( vtkFloatingPointType theValue )
66 {
67   if ( VISU::CheckIsSameValue( theValue, myShrinkFactor ) )
68     return;
69
70   myShrinkFactor = theValue;
71   this->Modified();
72 }
73
74
75 //----------------------------------------------------------------------------
76 vtkFloatingPointType VISU_ElnoDisassembleFilter::GetShrinkFactor()
77 {
78   return myShrinkFactor;
79 }
80
81
82 //----------------------------------------------------------------------------
83 namespace
84 {
85   //----------------------------------------------------------------------------
86   template < int points_type, int elno_type >
87   struct TExecute2
88   {
89     vtkUnstructuredGrid *myInput;
90     vtkUnstructuredGrid *myOutput;
91     vtkDataArray *myElnoDataArray;
92     vtkDataArray *myElnoDataMapper;
93     vtkFloatingPointType myShrinkFactor;
94
95     typedef typename VISU::TL::TEnum2VTKArrayType< points_type >::TResult TPointsDataArray;
96     typedef typename VISU::TL::TEnum2VTKBasicType< points_type >::TResult TPointsDataType;
97
98     typedef typename VISU::TL::TEnum2VTKArrayType< elno_type >::TResult TElnoDataArray;
99     typedef typename VISU::TL::TEnum2VTKBasicType< elno_type >::TResult TElnoDataType;
100
101     VISU::TGetElnoNodeData< elno_type > myGetElnoNodeData;
102     vtkCellArray *myConnectivity;
103     vtkPointData *myInputPointData;
104     vtkPointData *myOutputPointData;
105     TPointsDataArray *myInputPointsArray;
106     TPointsDataArray *myOutputPointsArray;
107     TElnoDataArray* myElnoFullDataArray;
108     TElnoDataArray* myElnoPartialDataArray;
109     TPointsDataArray *myElnoPointCoords;
110     vtkIntArray* myInputPointsMapper;
111     vtkIntArray* myOutputPointsMapper;
112
113     //----------------------------------------------------------------------------
114     TExecute2( vtkUnstructuredGrid *theInput, 
115                vtkUnstructuredGrid *theOutput, 
116                vtkDataArray *theElnoDataArray,
117                vtkDataArray *theElnoDataMapper,
118                vtkFloatingPointType theShrinkFactor )
119       : myGetElnoNodeData( theElnoDataArray, theElnoDataMapper )
120       , myInput( theInput )
121       , myOutput( theOutput )
122       , myElnoDataArray( theElnoDataArray )
123       , myElnoDataMapper( theElnoDataMapper )
124       , myShrinkFactor( theShrinkFactor )
125     {
126       myConnectivity = vtkCellArray::New();
127       myConnectivity->DeepCopy( theInput->GetCells() );
128     
129       vtkPoints *anInputPoints = theInput->GetPoints();
130       vtkPoints *aPoints = anInputPoints->New( anInputPoints->GetDataType() );
131       vtkIdType aNbCells = myConnectivity->GetNumberOfCells();
132       vtkIdType aNbPoints = myConnectivity->GetNumberOfConnectivityEntries() - aNbCells;
133       aPoints->Allocate( aNbPoints );
134     
135       myInputPointsArray = TPointsDataArray::SafeDownCast( anInputPoints->GetData() );
136       myOutputPointsArray = TPointsDataArray::SafeDownCast( aPoints->GetData() );
137     
138       myInputPointData = theInput->GetPointData();
139       myOutputPointData = theOutput->GetPointData();
140       myOutputPointData->Allocate( aNbPoints );
141     
142       vtkCellData *anInputCellData = theInput->GetCellData();
143
144       // To create a new copy of initial data for output
145       myElnoFullDataArray = TElnoDataArray::New();
146       myElnoFullDataArray->SetName( "VISU_FIELD" );
147       myElnoFullDataArray->SetNumberOfComponents( myGetElnoNodeData.getNbComp() );
148       myElnoFullDataArray->SetNumberOfTuples( aNbPoints );
149
150       // To create a new copy of partial initial data for output
151       myElnoPartialDataArray = TElnoDataArray::New();
152       // This partial data can be represented as in terms of vectors as scalars
153       if ( anInputCellData->GetVectors() != NULL ) 
154         myElnoPartialDataArray->SetNumberOfComponents( 3 );
155       else
156         myElnoPartialDataArray->SetNumberOfComponents( 1 );
157       myElnoPartialDataArray->SetNumberOfTuples( aNbPoints );
158
159       myElnoPointCoords = TPointsDataArray::New();
160       myElnoPointCoords->SetName( "ELNO_POINT_COORDS" );
161       myElnoPointCoords->SetNumberOfComponents( 3 );
162       myElnoPointCoords->SetNumberOfTuples( aNbPoints );
163
164       vtkDataArray* anArray = myInputPointData->GetArray( "VISU_POINTS_MAPPER" );
165       myInputPointsMapper = vtkIntArray::SafeDownCast( anArray );
166       
167       myOutputPointsMapper = vtkIntArray::New();
168       myOutputPointsMapper->SetName( myInputPointsMapper->GetName() );
169       myOutputPointsMapper->SetNumberOfComponents( myInputPointsMapper->GetNumberOfComponents() );
170       myOutputPointsMapper->SetNumberOfTuples( aNbPoints );
171
172       if ( theShrinkFactor > 0.0 )
173         this->ShrinkExecute();
174       else
175         this->SimpleExecute();
176
177       theOutput->SetPoints( aPoints );
178       
179       theOutput->SetCells( theInput->GetCellTypesArray(), 
180                            theInput->GetCellLocationsArray(),
181                            myConnectivity );
182
183       myConnectivity->Delete();
184       
185       vtkCellData *anOutputCellData = theOutput->GetCellData();
186       anOutputCellData->PassData( anInputCellData );
187       
188       anOutputCellData->RemoveArray( "ELNO_COMPONENT_MAPPER" );
189       anOutputCellData->RemoveArray( "ELNO_FIELD" );
190       anOutputCellData->RemoveArray( "VISU_FIELD" );
191       anOutputCellData->SetVectors( NULL );
192       
193       //anOutputPointData->PassData( anInputPointData );
194       
195       myOutputPointData->AddArray( myElnoFullDataArray );
196       myElnoFullDataArray->Delete();
197       
198       if ( anInputCellData->GetVectors() != NULL ) 
199         myOutputPointData->SetVectors( myElnoPartialDataArray );
200       else
201         myOutputPointData->SetScalars( myElnoPartialDataArray );
202       myElnoPartialDataArray->Delete();
203       
204       myOutputPointData->AddArray( myElnoPointCoords );
205       myElnoPointCoords->Delete();
206       
207       myOutputPointData->AddArray( myOutputPointsMapper );
208       myOutputPointsMapper->Delete();
209     }
210
211     //----------------------------------------------------------------------------
212     void SimpleExecute()
213     {
214       // To reserve a temproary value holder
215       vtkIdType aNbComp = std::max( 3, myGetElnoNodeData.getNbComp() );
216       std::vector< TElnoDataType > anElnoDataValues( aNbComp ); 
217
218       std::vector< int > anPointsMapperValues( myInputPointsMapper->GetNumberOfComponents() ); 
219    
220       myConnectivity->InitTraversal();
221       vtkIdType aNbPts = 0, *aPts = 0;
222       for ( vtkIdType aCellId = 0; myConnectivity->GetNextCell( aNbPts, aPts ); aCellId++ ) {
223         for ( vtkIdType aPntId = 0; aPntId < aNbPts; aPntId++ ) {
224           TPointsDataType aCoords[ 3 ];
225           vtkIdType aCurrentPntId = aPts[ aPntId ];
226           myInputPointsArray->GetTupleValue( aCurrentPntId, aCoords );
227           
228           aPts[ aPntId ] = myOutputPointsArray->InsertNextTupleValue( aCoords );
229           vtkIdType aNewPntId = aPts[ aPntId ];
230           
231           myElnoPointCoords->SetTupleValue( aNewPntId, aCoords );
232           
233           myOutputPointData->CopyData( myInputPointData, aCurrentPntId, aNewPntId );
234           
235           TElnoDataType* anElnoData = myGetElnoNodeData( aCellId, aPntId );
236           myElnoFullDataArray->SetTupleValue( aNewPntId,  anElnoData );
237           
238           myElnoFullDataArray->GetTupleValue( aNewPntId, &anElnoDataValues[ 0 ] );
239           myElnoPartialDataArray->SetTupleValue( aNewPntId, &anElnoDataValues[ 0 ] );
240
241           myInputPointsMapper->GetTupleValue( aCurrentPntId, &anPointsMapperValues[ 0 ] );
242           myOutputPointsMapper->SetTupleValue( aNewPntId, &anPointsMapperValues[ 0 ] );
243         }
244       }
245     }
246
247     //----------------------------------------------------------------------------
248     void ShrinkExecute()
249     {
250       // To reserve a temproary value holder
251       vtkIdType aNbComp = std::max( 3, myGetElnoNodeData.getNbComp() );
252       std::vector< TElnoDataType > anElnoDataValues( aNbComp ); 
253       
254       std::vector< int > anPointsMapperValues( myInputPointsMapper->GetNumberOfComponents() ); 
255    
256       myConnectivity->InitTraversal();
257       vtkIdType aNbPts = 0, *aPts = 0;
258       for ( vtkIdType aCellId = 0; myConnectivity->GetNextCell( aNbPts, aPts ); aCellId++ ) {
259         
260         TPointsDataType aCenter[ 3 ] = { TPointsDataType(), TPointsDataType(), TPointsDataType() };
261         
262         for ( vtkIdType aPntId = 0; aPntId < aNbPts; aPntId++ ) {
263           TPointsDataType aCoords[ 3 ];
264           myInputPointsArray->GetTupleValue( aPts[ aPntId ], aCoords );
265           
266           aCenter[ 0 ] += aCoords[ 0 ];
267           aCenter[ 1 ] += aCoords[ 1 ];
268           aCenter[ 2 ] += aCoords[ 2 ];
269         }
270         
271         aCenter[ 0 ] /= aNbPts;
272         aCenter[ 1 ] /= aNbPts;
273         aCenter[ 2 ] /= aNbPts;
274         
275         for ( vtkIdType aPntId = 0; aPntId < aNbPts; aPntId++ ) {
276           TPointsDataType aCoords[ 3 ];
277           vtkIdType aCurrentPntId = aPts[ aPntId ];
278           myInputPointsArray->GetTupleValue( aCurrentPntId, aCoords );
279           
280           TPointsDataType aNewCoords[ 3 ];
281           
282           aNewCoords[ 0 ] = aCenter[ 0 ]  + 
283             TPointsDataType( myShrinkFactor * ( aCoords[ 0 ] - aCenter[ 0 ] ) );
284           aNewCoords[ 1 ] = aCenter[ 1 ]  + 
285             TPointsDataType( myShrinkFactor * ( aCoords[ 1 ] - aCenter[ 1 ] ) );
286           aNewCoords[ 2 ] = aCenter[ 2 ]  + 
287             TPointsDataType( myShrinkFactor * ( aCoords[ 2 ] - aCenter[ 2 ] ) );
288           
289           aPts[ aPntId ] = myOutputPointsArray->InsertNextTupleValue( aNewCoords );
290           vtkIdType aNewPntId = aPts[ aPntId ];
291           
292           myElnoPointCoords->SetTupleValue( aNewPntId, aCoords );
293           
294           myOutputPointData->CopyData( myInputPointData, aCurrentPntId, aNewPntId );
295           
296           TElnoDataType* anElnoData = myGetElnoNodeData( aCellId, aPntId );
297           myElnoFullDataArray->SetTupleValue( aNewPntId, anElnoData );
298           
299           myElnoFullDataArray->GetTupleValue( aNewPntId, &anElnoDataValues[ 0 ] );
300           myElnoPartialDataArray->SetTupleValue( aNewPntId, &anElnoDataValues[ 0 ] );
301
302           myInputPointsMapper->GetTupleValue( aCurrentPntId, &anPointsMapperValues[ 0 ] );
303           myOutputPointsMapper->SetTupleValue( aNewPntId, &anPointsMapperValues[ 0 ] );
304         }
305       }
306     }
307   };
308
309
310   //----------------------------------------------------------------------------
311   template < int points_type, int elno_type >
312   int Execute2( vtkUnstructuredGrid *theInput, 
313                 vtkUnstructuredGrid *theOutput, 
314                 vtkDataArray *theElnoDataArray,
315                 vtkDataArray *theElnoDataMapper,
316                 vtkFloatingPointType theShrinkFactor )
317   {
318     TExecute2< points_type, elno_type >( theInput, 
319                                          theOutput, 
320                                          theElnoDataArray, 
321                                          theElnoDataMapper, 
322                                          theShrinkFactor );
323         
324     return 1;
325   }
326           
327
328   //----------------------------------------------------------------------------
329   template < int points_type >
330   int Execute( vtkUnstructuredGrid *theInput, 
331                vtkUnstructuredGrid *theOutput, 
332                vtkDataArray *theElnoDataArray,
333                vtkDataArray *theElnoDataMapper,
334                vtkFloatingPointType theShrinkFactor )
335   {
336     switch( theElnoDataArray->GetDataType() ){
337     case VTK_DOUBLE:
338       return Execute2< points_type, VTK_DOUBLE >
339         ( theInput, theOutput, theElnoDataArray, theElnoDataMapper, theShrinkFactor );
340     case VTK_FLOAT:
341       return Execute2< points_type, VTK_FLOAT >
342         ( theInput, theOutput, theElnoDataArray, theElnoDataMapper, theShrinkFactor );
343     case VTK_INT:
344       return Execute2< points_type, VTK_INT >
345         ( theInput, theOutput, theElnoDataArray, theElnoDataMapper, theShrinkFactor );
346     case VTK_LONG:
347       return Execute2< points_type, VTK_LONG >
348         ( theInput, theOutput, theElnoDataArray, theElnoDataMapper, theShrinkFactor );
349     default:
350       break;
351     }
352     
353     return 0;
354   } 
355
356
357   //----------------------------------------------------------------------------
358 }
359
360
361 //----------------------------------------------------------------------------
362 int VISU_ElnoDisassembleFilter::RequestData( vtkInformation *vtkNotUsed(request),
363                                              vtkInformationVector **inputVector,
364                                              vtkInformationVector *outputVector )
365 {
366   // get the info objects
367   vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
368   vtkInformation *outInfo = outputVector->GetInformationObject(0);
369
370   // get the input and ouptut
371   vtkUnstructuredGrid *anInput =
372     vtkUnstructuredGrid::SafeDownCast( inInfo->Get( vtkDataObject::DATA_OBJECT() ) );
373   vtkUnstructuredGrid *anOutput = 
374     vtkUnstructuredGrid::SafeDownCast( outInfo->Get( vtkDataObject::DATA_OBJECT() ) );
375
376   vtkDataArray *anElnoDataArray = this->GetInputArrayToProcess( 0, inputVector );
377   vtkDataArray *anElnoDataMapper = this->GetInputArrayToProcess( 1, inputVector );
378
379   if ( !anElnoDataArray ) {
380     anOutput->ShallowCopy( anInput );
381     return 1;
382   }
383
384   vtkPoints *aPoints = anInput->GetPoints();
385   switch( aPoints->GetDataType() ){
386   case VTK_DOUBLE:
387     return ::Execute< VTK_DOUBLE >( anInput, anOutput, anElnoDataArray, anElnoDataMapper, myShrinkFactor );
388   case VTK_FLOAT:
389     return ::Execute< VTK_FLOAT >( anInput, anOutput, anElnoDataArray, anElnoDataMapper, myShrinkFactor );
390   case VTK_INT:
391     return ::Execute< VTK_INT >( anInput, anOutput, anElnoDataArray, anElnoDataMapper, myShrinkFactor );
392   case VTK_LONG:
393     return ::Execute< VTK_LONG >( anInput, anOutput, anElnoDataArray, anElnoDataMapper, myShrinkFactor );
394   default:
395     break;
396   }  
397   
398   return 0;
399 }
400
401
402 //----------------------------------------------------------------------------