Salome HOME
de46b628a0d29019d7283e50cc30420e7d50e50c
[modules/visu.git] / src / PIPELINE / VISU_ElnoAssembleFilter.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_ElnoAssembleFilter.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_ElnoAssembleFilter );
38
39
40 //----------------------------------------------------------------------------
41 VISU_ElnoAssembleFilter::VISU_ElnoAssembleFilter()
42 {
43   this->SetInputArrayToProcess( 0, // idx
44                                 0, // port
45                                 0, // connection
46                                 vtkDataObject::FIELD_ASSOCIATION_POINTS, // field association
47                                 "ELNO_POINT_COORDS" ); // name
48
49   this->myIsRestorePoints = false;
50 }
51
52
53 //----------------------------------------------------------------------------
54 VISU_ElnoAssembleFilter::~VISU_ElnoAssembleFilter()
55 {}
56
57
58 //----------------------------------------------------------------------------
59 void VISU_ElnoAssembleFilter::SetElnoAssembleState( bool theIsRestorePoints )
60 {
61   if ( myIsRestorePoints == theIsRestorePoints )
62     return;
63     
64   myIsRestorePoints = theIsRestorePoints;
65   this->Modified();
66 }
67
68 //----------------------------------------------------------------------------
69 namespace
70 {
71   //----------------------------------------------------------------------------
72   template < int points_type, int elno_type >
73   int Execute2( vtkPointSet *theInput, 
74                 vtkPointSet *theOutput,
75                 vtkDataArray *theElnoPointCoords )
76   {
77     theOutput->CopyStructure( theInput );
78     
79     vtkCellData *aCellData = theOutput->GetCellData();
80     aCellData->PassData( theInput->GetCellData() );
81
82     vtkPointData *aPointData = theOutput->GetPointData();
83     aPointData->PassData( theInput->GetPointData() );
84
85     vtkPoints *anInputPoints = theInput->GetPoints();
86     vtkPoints *aPoints = anInputPoints->New( elno_type );
87     vtkIdType aNbPoints = theInput->GetNumberOfPoints();
88     aPoints->SetNumberOfPoints( aNbPoints );
89     
90     typedef typename VISU::TL::TEnum2VTKArrayType< elno_type >::TResult TPointsDataArray;
91     typedef typename VISU::TL::TEnum2VTKBasicType< elno_type >::TResult TPointsDataType;
92     TPointsDataArray* anOutputPointsArray = TPointsDataArray::SafeDownCast( aPoints->GetData() );
93
94     TPointsDataArray* anElnoPointCoords = TPointsDataArray::SafeDownCast( theElnoPointCoords );
95     
96     for ( vtkIdType aPointId = 0; aPointId < aNbPoints; aPointId++ ) {
97       TPointsDataType aCoords[ 3 ];
98       anElnoPointCoords->GetTupleValue( aPointId, aCoords );
99       anOutputPointsArray->SetTupleValue( aPointId, aCoords );
100     }
101     
102     theOutput->SetPoints( aPoints );
103
104     return 1;
105   } 
106
107
108   //----------------------------------------------------------------------------
109   template < int points_type >
110   int Execute( vtkPointSet *theInput, 
111                vtkPointSet *theOutput,
112                vtkDataArray *theElnoPointCoords )
113   {
114     switch( theElnoPointCoords->GetDataType() ){
115     case VTK_DOUBLE:
116       return Execute2< points_type, VTK_DOUBLE >( theInput, theOutput, theElnoPointCoords );
117     case VTK_FLOAT:
118       return Execute2< points_type, VTK_FLOAT >( theInput, theOutput, theElnoPointCoords );
119     case VTK_INT:
120       return Execute2< points_type, VTK_INT >( theInput, theOutput, theElnoPointCoords );
121     case VTK_LONG:
122       return Execute2< points_type, VTK_LONG >( theInput, theOutput, theElnoPointCoords );
123     default:
124       break;
125     }
126     
127     return 0;
128   } 
129
130
131   //----------------------------------------------------------------------------
132 }
133
134
135 //----------------------------------------------------------------------------
136 int VISU_ElnoAssembleFilter::RequestData( vtkInformation *vtkNotUsed(request),
137                                           vtkInformationVector **inputVector,
138                                           vtkInformationVector *outputVector )
139 {
140   // get the info objects
141   vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
142   vtkInformation *outInfo = outputVector->GetInformationObject(0);
143
144   // get the input and ouptut
145   vtkPointSet *anInput = vtkPointSet::SafeDownCast( inInfo->Get( vtkDataObject::DATA_OBJECT() ) );
146   vtkPointSet *anOutput = vtkPointSet::SafeDownCast( outInfo->Get( vtkDataObject::DATA_OBJECT() ) );
147
148   vtkDataArray *anElnoPointCoords = this->GetInputArrayToProcess( 0, inputVector );
149
150   if ( !myIsRestorePoints || !anElnoPointCoords ) {
151     anOutput->ShallowCopy( anInput );
152     return 1;
153   }
154
155   vtkPoints *aPoints = anInput->GetPoints();
156   switch( aPoints->GetDataType() ){
157   case VTK_DOUBLE:
158     return ::Execute< VTK_DOUBLE >( anInput, anOutput, anElnoPointCoords );
159   case VTK_FLOAT:
160     return ::Execute< VTK_FLOAT >( anInput, anOutput, anElnoPointCoords );
161   case VTK_INT:
162     return ::Execute< VTK_INT >( anInput, anOutput, anElnoPointCoords );
163   case VTK_LONG:
164     return ::Execute< VTK_LONG >( anInput, anOutput, anElnoPointCoords );
165   default:
166     break;
167   }  
168   
169   return 0;
170 }
171
172
173 //----------------------------------------------------------------------------