Salome HOME
ea9a5c0ac790501af74d3530ba1d82179ceb8a45
[modules/visu.git] / src / PIPELINE / VISU_ElnoAssembleFilter.cxx
1 // Copyright (C) 2007-2011  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "VISU_ElnoAssembleFilter.hxx"
21 #include "VISU_PipeLineUtils.hxx"
22 #include "VISU_ElnoMeshValue.hxx"
23
24 #include <vtkCellData.h>
25 #include <vtkInformation.h>
26 #include <vtkInformationVector.h>
27 #include <vtkObjectFactory.h>
28 #include <vtkPointData.h>
29 #include <vtkUnstructuredGrid.h>
30 #include <vtkPoints.h>
31 #include <vtkCellArray.h>
32
33
34 //----------------------------------------------------------------------------
35 vtkStandardNewMacro( VISU_ElnoAssembleFilter );
36
37
38 //----------------------------------------------------------------------------
39 VISU_ElnoAssembleFilter::VISU_ElnoAssembleFilter()
40 {
41   this->SetInputArrayToProcess( 0, // idx
42                                 0, // port
43                                 0, // connection
44                                 vtkDataObject::FIELD_ASSOCIATION_POINTS, // field association
45                                 "ELNO_POINT_COORDS" ); // name
46
47   this->myIsRestorePoints = false;
48 }
49
50
51 //----------------------------------------------------------------------------
52 VISU_ElnoAssembleFilter::~VISU_ElnoAssembleFilter()
53 {}
54
55
56 //----------------------------------------------------------------------------
57 void VISU_ElnoAssembleFilter::SetElnoAssembleState( bool theIsRestorePoints )
58 {
59   if ( myIsRestorePoints == theIsRestorePoints )
60     return;
61     
62   myIsRestorePoints = theIsRestorePoints;
63   this->Modified();
64 }
65
66 //----------------------------------------------------------------------------
67 namespace
68 {
69   //----------------------------------------------------------------------------
70   template < int points_type, int elno_type >
71   int Execute2( vtkPointSet *theInput, 
72                 vtkPointSet *theOutput,
73                 vtkDataArray *theElnoPointCoords )
74   {
75     theOutput->CopyStructure( theInput );
76     
77     vtkCellData *aCellData = theOutput->GetCellData();
78     aCellData->PassData( theInput->GetCellData() );
79
80     vtkPointData *aPointData = theOutput->GetPointData();
81     aPointData->PassData( theInput->GetPointData() );
82
83     vtkPoints *anInputPoints = theInput->GetPoints();
84     vtkPoints *aPoints = anInputPoints->New( elno_type );
85     vtkIdType aNbPoints = theInput->GetNumberOfPoints();
86     aPoints->SetNumberOfPoints( aNbPoints );
87     
88     typedef typename VISU::TL::TEnum2VTKArrayType< elno_type >::TResult TPointsDataArray;
89     typedef typename VISU::TL::TEnum2VTKBasicType< elno_type >::TResult TPointsDataType;
90     TPointsDataArray* anOutputPointsArray = TPointsDataArray::SafeDownCast( aPoints->GetData() );
91
92     TPointsDataArray* anElnoPointCoords = TPointsDataArray::SafeDownCast( theElnoPointCoords );
93     
94     for ( vtkIdType aPointId = 0; aPointId < aNbPoints; aPointId++ ) {
95       TPointsDataType aCoords[ 3 ];
96       anElnoPointCoords->GetTupleValue( aPointId, aCoords );
97       anOutputPointsArray->SetTupleValue( aPointId, aCoords );
98     }
99     
100     theOutput->SetPoints( aPoints );
101
102     return 1;
103   } 
104
105
106   //----------------------------------------------------------------------------
107   template < int points_type >
108   int Execute( vtkPointSet *theInput, 
109                vtkPointSet *theOutput,
110                vtkDataArray *theElnoPointCoords )
111   {
112     switch( theElnoPointCoords->GetDataType() ){
113     case VTK_DOUBLE:
114       return Execute2< points_type, VTK_DOUBLE >( theInput, theOutput, theElnoPointCoords );
115     case VTK_FLOAT:
116       return Execute2< points_type, VTK_FLOAT >( theInput, theOutput, theElnoPointCoords );
117     case VTK_INT:
118       return Execute2< points_type, VTK_INT >( theInput, theOutput, theElnoPointCoords );
119     case VTK_LONG:
120       return Execute2< points_type, VTK_LONG >( theInput, theOutput, theElnoPointCoords );
121     default:
122       break;
123     }
124     
125     return 0;
126   } 
127
128
129   //----------------------------------------------------------------------------
130 }
131
132
133 //----------------------------------------------------------------------------
134 int VISU_ElnoAssembleFilter::RequestData( vtkInformation *vtkNotUsed(request),
135                                           vtkInformationVector **inputVector,
136                                           vtkInformationVector *outputVector )
137 {
138   // get the info objects
139   vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
140   vtkInformation *outInfo = outputVector->GetInformationObject(0);
141
142   // get the input and ouptut
143   vtkPointSet *anInput = vtkPointSet::SafeDownCast( inInfo->Get( vtkDataObject::DATA_OBJECT() ) );
144   vtkPointSet *anOutput = vtkPointSet::SafeDownCast( outInfo->Get( vtkDataObject::DATA_OBJECT() ) );
145
146   vtkDataArray *anElnoPointCoords = this->GetInputArrayToProcess( 0, inputVector );
147
148   if ( !myIsRestorePoints || !anElnoPointCoords ) {
149     anOutput->ShallowCopy( anInput );
150     return 1;
151   }
152
153   vtkPoints *aPoints = anInput->GetPoints();
154   switch( aPoints->GetDataType() ){
155   case VTK_DOUBLE:
156     return ::Execute< VTK_DOUBLE >( anInput, anOutput, anElnoPointCoords );
157   case VTK_FLOAT:
158     return ::Execute< VTK_FLOAT >( anInput, anOutput, anElnoPointCoords );
159   case VTK_INT:
160     return ::Execute< VTK_INT >( anInput, anOutput, anElnoPointCoords );
161   case VTK_LONG:
162     return ::Execute< VTK_LONG >( anInput, anOutput, anElnoPointCoords );
163   default:
164     break;
165   }  
166   
167   return 0;
168 }
169
170
171 //----------------------------------------------------------------------------