Salome HOME
105e7f5f4f45238cae553fd17a66e586ad722e23
[modules/visu.git] / src / PIPELINE / VISU_Extractor.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 //  VISU OBJECT : interactive object for VISU entities implementation
23 //  File   : VISU_Extractor.cxx
24 //  Module : VISU
25
26 #include "VISU_Extractor.hxx"
27 #include "VISU_PipeLineUtils.hxx"
28 #include "VISU_ConvertorUtils.hxx"
29 #include "VISU_MeshValue.hxx"
30
31 #include <sstream>
32
33 #include <vtkObjectFactory.h>
34 #include <vtkUnstructuredGrid.h>
35 #include <vtkPointData.h>
36 #include <vtkCellData.h>
37 #include <vector>
38 #include <vtkInformation.h>
39 #include <vtkInformationVector.h>
40
41
42 //----------------------------------------------------------------------------
43 vtkStandardNewMacro(VISU_Extractor);
44
45 //----------------------------------------------------------------------------
46 VISU_Extractor
47 ::VISU_Extractor():
48   myScalarMode(1),
49   myGaussMetric(VISU::AVERAGE_METRIC)
50 {
51 }
52
53
54 //----------------------------------------------------------------------------
55 VISU_Extractor
56 ::~VISU_Extractor()
57 {}
58
59
60 //----------------------------------------------------------------------------
61 void
62 VISU_Extractor
63 ::SetScalarMode(int theScalarMode)
64 {
65   if(myScalarMode != theScalarMode){
66     myScalarMode = theScalarMode;
67     Modified();
68   }
69 }
70
71 //----------------------------------------------------------------------------
72 int
73 VISU_Extractor
74 ::GetScalarMode()
75 {
76   return myScalarMode;
77 }
78
79
80 //----------------------------------------------------------------------------
81 void
82 VISU_Extractor
83 ::SetGaussMetric(VISU::TGaussMetric theGaussMetric)
84 {
85   if(myGaussMetric != theGaussMetric){
86     myGaussMetric = theGaussMetric;
87     Modified();
88   }
89 }
90
91 //----------------------------------------------------------------------------
92 VISU::TGaussMetric
93 VISU_Extractor
94 ::GetGaussMetric()
95 {
96   return myGaussMetric;
97 }
98
99
100 //----------------------------------------------------------------------------
101 template<typename TValueType>
102 void
103 Module2Scalars(vtkDataArray *theInputDataArray,
104                TValueType* theOutputPtr,
105                vtkIdType theNbOfTuples,
106                VISU::TGaussMetric theGaussMetric)
107 {
108   vtkIdType aNbComp = theInputDataArray->GetNumberOfComponents();
109   if (aNbComp != 3) // Min, Max, Avg
110     return;
111   std::vector<vtkFloatingPointType> anArray (3);
112   for (vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++) {
113     theInputDataArray->GetTuple(aTupleId, &anArray[0]);
114     switch (theGaussMetric) {
115       case VISU::MINIMUM_METRIC: *theOutputPtr = TValueType(anArray[0]); break;
116       case VISU::MAXIMUM_METRIC: *theOutputPtr = TValueType(anArray[1]); break;
117       case VISU::AVERAGE_METRIC: *theOutputPtr = TValueType(anArray[2]); break;
118     }
119     theOutputPtr++;
120   }
121 }
122
123
124 //----------------------------------------------------------------------------
125 template<typename TValueType>
126 void
127 Component2Scalars(vtkDataArray *theInputDataArray,
128                   TValueType* theInputPtr,
129                   TValueType* theOutputPtr,
130                   vtkIdType theNbOfTuples,
131                   vtkIdType theComponentId)
132 {
133   vtkIdType aNbComp = theInputDataArray->GetNumberOfComponents();
134   for(vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++){
135     *theOutputPtr = *(theInputPtr + theComponentId);
136     theInputPtr += aNbComp;
137     theOutputPtr++;
138   }
139 }
140
141
142 //----------------------------------------------------------------------------
143 template<typename TDataSetAttributesType> void
144 ExecuteScalars(vtkIdType theNbOfTuples,
145                vtkIdType theScalarMode,
146                VISU::TGaussMetric theGaussMetric,
147                TDataSetAttributesType* theInputData,
148                TDataSetAttributesType* theOutputData)
149 {
150   if (theNbOfTuples < 1)
151     return;
152
153   vtkDataArray* aFieldArray = NULL;
154   switch (theGaussMetric) {
155     case VISU::AVERAGE_METRIC: aFieldArray = theInputData->GetArray("VISU_FIELD"); break;
156     case VISU::MINIMUM_METRIC: aFieldArray = theInputData->GetArray("VISU_FIELD_GAUSS_MIN"); break;
157     case VISU::MAXIMUM_METRIC: aFieldArray = theInputData->GetArray("VISU_FIELD_GAUSS_MAX"); break;
158   }
159   if( !aFieldArray )
160     return;
161
162   vtkIdType anInputDataType = aFieldArray->GetDataType();
163   vtkDataArray *anOutputScalars = vtkDataArray::CreateDataArray(anInputDataType);
164   anOutputScalars->SetNumberOfComponents(1);
165   anOutputScalars->SetNumberOfTuples(theNbOfTuples);
166
167   void *anInputPtr = aFieldArray->GetVoidPointer(0);
168   void *anOutputPtr = anOutputScalars->GetVoidPointer(0);
169
170   if (theScalarMode == 0) {
171     aFieldArray = theInputData->GetArray("VISU_FIELD_GAUSS_MOD");
172     if (!aFieldArray)
173       return;
174     switch (anInputDataType) {
175       vtkTemplateMacro4(Module2Scalars,
176                         aFieldArray,
177                         (VTK_TT *)(anOutputPtr),
178                         theNbOfTuples,
179                         theGaussMetric);
180     default:
181       break;
182     }
183   } else {
184     switch (anInputDataType) {
185       vtkTemplateMacro5(Component2Scalars,
186                         aFieldArray,
187                         (VTK_TT *)(anInputPtr),
188                         (VTK_TT *)(anOutputPtr),
189                         theNbOfTuples,
190                         theScalarMode - 1);
191     default:
192       break;
193     }
194   }
195
196   theOutputData->SetScalars(anOutputScalars);
197   anOutputScalars->Delete();
198 }
199
200
201 //---------------------------------------------------------------
202 int
203 VISU_Extractor
204 ::RequestData(vtkInformation *theRequest,
205               vtkInformationVector **theInputVector,
206               vtkInformationVector *theOutputVector)
207 {
208   vtkDataSet *anInput = VISU::GetInput( theInputVector, 0 );
209   vtkDataSet *anOutput = VISU::GetOutput( theOutputVector );
210
211   anOutput->CopyStructure( anInput );
212
213   vtkPointData *anInputPointData = anInput->GetPointData();
214   vtkPointData *anOutputPointData = anOutput->GetPointData();
215   anOutputPointData->PassData( anInputPointData );
216   if ( VISU::IsDataOnPoints( anInput ) ) {
217     int aNbElems = anInput->GetNumberOfPoints();
218     if ( anInputPointData->GetAttribute( vtkDataSetAttributes::VECTORS ) )
219       ExecuteScalars( aNbElems, myScalarMode, myGaussMetric, anInputPointData, anOutputPointData );
220   }
221
222   vtkCellData *anInputCellData = anInput->GetCellData();
223   vtkCellData *anOutputCellData = anOutput->GetCellData();
224   anOutputCellData->PassData( anInputCellData );
225   if ( VISU::IsDataOnCells( anInput ) ) {
226     int aNbElems = anInput->GetNumberOfCells();
227     if ( anInputCellData->GetAttribute( vtkDataSetAttributes::VECTORS ) )
228       ExecuteScalars( aNbElems, myScalarMode, myGaussMetric, anInputCellData, anOutputCellData );
229   }
230
231   return 1;
232 }