]> SALOME platform Git repositories - modules/visu.git/blob - src/PIPELINE/VISU_ColoredPL.cxx
Salome HOME
windows port for isnan
[modules/visu.git] / src / PIPELINE / VISU_ColoredPL.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_ColoredPL.cxx
24 // Author:  Alexey PETROV
25 // Module : VISU
26 //
27 #include "VISU_ColoredPL.hxx"
28 #include "VISU_Extractor.hxx"
29 #include "VISU_FieldTransform.hxx"
30 #include "VISU_LookupTable.hxx"
31 #include "VISU_MapperHolder.hxx"
32
33 #include "VISU_PipeLineUtils.hxx"
34
35 #include <vtkThreshold.h>
36 #include <vtkPassThroughFilter.h>
37 #include <vtkDoubleArray.h>
38
39 #ifdef WNT
40 #include <float.h>
41 #define isnan _isnan
42 #endif
43
44 //----------------------------------------------------------------------------
45 VISU_ColoredPL
46 ::VISU_ColoredPL():
47   myMapperTable( VISU_LookupTable::New() ),
48   myBarTable( VISU_LookupTable::New() ),
49   myExtractor( VISU_Extractor::New() ),
50   myFieldTransform( VISU_FieldTransform::New() ),
51   myThreshold ( vtkThreshold::New() ),
52   myPassFilter( vtkPassThroughFilter::New() ),
53   myDistribution( vtkDoubleArray::New() )
54 {
55   myMapperTable->Delete();
56   myMapperTable->SetScale(VTK_SCALE_LINEAR);
57   myMapperTable->SetHueRange(0.667, 0.0);
58
59   myBarTable->Delete();
60   myBarTable->SetScale(VTK_SCALE_LINEAR);
61   myBarTable->SetHueRange(0.667, 0.0);
62
63   myExtractor->Delete();
64
65   myFieldTransform->Delete();
66
67   myThreshold->AllScalarsOn(); 
68   myThreshold->Delete();
69   myPassFilter->Delete();
70   myDistribution->Delete();
71 }
72
73
74 //----------------------------------------------------------------------------
75 VISU_ColoredPL
76 ::~VISU_ColoredPL()
77 {}
78
79
80 //----------------------------------------------------------------------------
81 unsigned long int 
82 VISU_ColoredPL
83 ::GetMTime()
84 {
85   unsigned long int aTime = Superclass::GetMTime();
86
87   aTime = std::max(aTime, myMapperTable->GetMTime());
88   aTime = std::max(aTime, myBarTable->GetMTime());
89   aTime = std::max(aTime, myExtractor->GetMTime());
90   aTime = std::max(aTime, myFieldTransform->GetMTime());
91   aTime = std::max(aTime, myThreshold->GetMTime());
92   aTime = std::max(aTime, myPassFilter->GetMTime());
93   aTime = std::max(aTime, myDistribution->GetMTime());
94
95   return aTime;
96 }
97
98
99 //----------------------------------------------------------------------------
100 void
101 VISU_ColoredPL
102 ::DoShallowCopy(VISU_PipeLine *thePipeLine,
103                 bool theIsCopyInput)
104 {
105   Superclass::DoShallowCopy(thePipeLine, theIsCopyInput);
106
107   if(VISU_ColoredPL *aPipeLine = dynamic_cast<VISU_ColoredPL*>(thePipeLine)){
108     if ( theIsCopyInput ) {
109       SetScalarRange( aPipeLine->GetScalarRange() );
110       if ( this->IsScalarFilterUsed() )
111         SetScalarFilterRange( aPipeLine->GetScalarFilterRange() );
112     }
113
114     SetScalarMode(aPipeLine->GetScalarMode());
115     SetNbColors(aPipeLine->GetNbColors());
116     SetScaling(aPipeLine->GetScaling());
117     SetMapScale(aPipeLine->GetMapScale());
118   }
119 }
120
121
122 //----------------------------------------------------------------------------
123 int
124 VISU_ColoredPL
125 ::GetScalarMode()
126 {
127   return myExtractor->GetScalarMode();
128 }
129
130
131 //----------------------------------------------------------------------------
132 void
133 VISU_ColoredPL
134 ::SetScalarMode(int theScalarMode,
135                 vtkDataSet *theInput,
136                 VISU_Extractor* theExtractor)
137 {
138   if(theInput){
139     if(VISU::IsDataOnPoints(theInput)){
140       vtkPointData *aPointData = theInput->GetPointData();
141       if(!aPointData->GetAttribute(vtkDataSetAttributes::VECTORS)) {
142         if(theScalarMode == 0){
143           return;
144         }
145       }
146     } else {
147       vtkCellData *aCellData = theInput->GetCellData();
148       if(!aCellData->GetAttribute(vtkDataSetAttributes::VECTORS)){
149         if(theScalarMode == 0){
150           return;
151         }
152       }
153     }
154   }
155
156   theExtractor->SetScalarMode(theScalarMode);
157 }
158
159 //----------------------------------------------------------------------------
160 void
161 VISU_ColoredPL
162 ::SetScalarMode(int theScalarMode)
163 {
164   SetScalarMode(theScalarMode, GetInput(), myExtractor);
165 }
166
167
168 //----------------------------------------------------------------------------
169 void
170 VISU_ColoredPL
171 ::SetScalarRange( vtkFloatingPointType theRange[2] )
172 {
173   if ( theRange[0] > theRange[1] ) 
174     return;
175   
176   if (VISU::CheckIsSameRange( GetScalarRange(), theRange) )
177     return;
178
179   myFieldTransform->SetScalarRange( theRange );
180   myBarTable->SetRange( theRange );
181 }
182
183
184 //----------------------------------------------------------------------------
185 void
186 VISU_ColoredPL
187 ::SetScalarFilterRange( vtkFloatingPointType theRange[2] )
188 {
189   vtkFloatingPointType aRange[ 2 ];
190   this->GetScalarFilterRange( aRange );
191
192   if ( VISU::CheckIsSameRange( aRange, theRange) )
193     return;
194
195   myThreshold->ThresholdBetween( theRange[0], theRange[1] );
196 }
197
198
199 //----------------------------------------------------------------------------
200 void
201 VISU_ColoredPL
202 ::GetScalarFilterRange( vtkFloatingPointType theRange[2] )
203 {
204   theRange[ 0 ] = myThreshold->GetLowerThreshold();
205   theRange[ 1 ] = myThreshold->GetUpperThreshold();
206 }
207
208
209 //----------------------------------------------------------------------------
210 vtkFloatingPointType*
211 VISU_ColoredPL
212 ::GetScalarFilterRange()
213 {
214   static vtkFloatingPointType aRange[ 2 ];
215
216   this->GetScalarFilterRange( aRange );
217
218   return aRange;
219 }
220
221
222 //----------------------------------------------------------------------------
223 void
224 VISU_ColoredPL
225 ::UseScalarFiltering( bool theUseScalarFilter )
226 {
227   if ( theUseScalarFilter ) {
228     // Include threshold filter between the transform and the pass filters. 
229     myPassFilter->SetInput( myThreshold->GetOutput() );
230   } else {
231     // Exclude threshold filter before the pass filter. 
232     myPassFilter->SetInput( myFieldTransform->GetOutput() );
233   }
234 }
235
236
237 //----------------------------------------------------------------------------
238 bool
239 VISU_ColoredPL
240 ::IsScalarFilterUsed()
241 {
242   return myThreshold->GetOutput() == myPassFilter->GetInput();
243 }
244
245
246 //----------------------------------------------------------------------------
247 vtkDoubleArray* 
248 VISU_ColoredPL
249 ::GetDistribution() 
250 {
251   unsigned long int aTime = this->GetMTime();
252   // If modified then update the distribution array
253   if (aTime > myDistribution->GetMTime()) {
254         // Set number of colors for the distribution
255     int nbColors = this->GetNbColors();
256         this->myDistribution->SetNumberOfValues(nbColors);
257         // Initialize numbers of colored cells with zero
258         this->myDistribution->FillComponent(0, 0);
259         // Create a lookup table to compute a color of a cell
260     VISU_LookupTable* lut = GetMapperTable();
261     vtkFloatingPointType aMapScale = lut->GetMapScale();
262     // Get scalar values from the input data to calculate their distribution within cells
263     vtkDataArray* dataArr;
264     // Dtermine where we have to take scalars from: cells data or points data. 
265     if(VISU::IsDataOnCells(this->GetOutput())) {
266         dataArr = this->GetOutput()->GetCellData()->GetScalars();
267     } else {
268         dataArr = this->GetOutput()->GetPointData()->GetScalars();
269     }
270     // If scalars data array is not defined then create an empty one to avoid exceptions
271     if (dataArr == NULL) {
272         dataArr = vtkDoubleArray::New();
273     }
274     
275     // Get range of scalars values
276 //    vtkFloatingPointType aRange[2];
277 //    dataArr->GetRange(aRange);
278
279     // Build the lookup table with the found range
280     // Get number of scalar values
281     int aNbVals = dataArr->GetNumberOfTuples();
282     if (aNbVals > 0) {
283       // Count the number of scalar values for each color in the input data
284       int idx = 0;
285       double cnt = 0;
286       // For each scalar value
287       for(vtkIdType aValId = 0; aValId < aNbVals; aValId++){
288         // Find the color index for this scalar value
289         idx = lut->GetIndex(*(dataArr->GetTuple(aValId)) * aMapScale);
290         // Increment the distribution value for this color index
291         cnt = this->myDistribution->GetValue(idx);
292         this->myDistribution->SetValue(idx, cnt + 1);
293       }
294       // Compute relative values when 1 is according to the total number of scalar values
295       for(vtkIdType aValId = 0; aValId < nbColors; aValId++){
296         cnt = this->myDistribution->GetValue(aValId);
297         this->myDistribution->SetValue(aValId, cnt / aNbVals);
298       }
299     }
300     this->myDistribution->Modified();
301         
302   }
303   
304   return myDistribution;
305 }
306 //----------------------------------------------------------------------------
307   // RKV : End
308
309 //----------------------------------------------------------------------------
310 vtkFloatingPointType* 
311 VISU_ColoredPL
312 ::GetScalarRange() 
313 {
314   return myFieldTransform->GetScalarRange();
315 }
316
317 //----------------------------------------------------------------------------
318 void
319 VISU_ColoredPL
320 ::SetScaling(int theScaling) 
321 {
322   if(GetScaling() == theScaling)
323     return;
324
325   myBarTable->SetScale(theScaling);
326
327   if(theScaling == VTK_SCALE_LOG10)
328     myFieldTransform->SetScalarTransform(&(VISU_FieldTransform::Log10));
329   else
330     myFieldTransform->SetScalarTransform(&(VISU_FieldTransform::Ident));
331 }
332
333 //----------------------------------------------------------------------------
334 int
335 VISU_ColoredPL
336 ::GetScaling() 
337 {
338   return myBarTable->GetScale();
339 }
340
341 //----------------------------------------------------------------------------
342 void
343 VISU_ColoredPL
344 ::SetNbColors(int theNbColors) 
345 {
346   myMapperTable->SetNumberOfColors(theNbColors);
347   myBarTable->SetNumberOfColors(theNbColors);
348 }
349
350 int
351 VISU_ColoredPL
352 ::GetNbColors() 
353 {
354   return myMapperTable->GetNumberOfColors();
355 }
356
357
358 //----------------------------------------------------------------------------
359 void
360 VISU_ColoredPL
361 ::Init()
362 {
363   SetScalarMode(0);
364
365   vtkFloatingPointType aRange[2];
366   GetSourceRange( aRange );
367
368   SetScalarRange( aRange );
369   SetScalarFilterRange( aRange );
370 }
371
372 //----------------------------------------------------------------------------
373 vtkPointSet* 
374 VISU_ColoredPL
375 ::GetClippedInput()
376 {
377   if(myPassFilter->GetInput())
378     myPassFilter->Update();
379   return myPassFilter->GetUnstructuredGridOutput();
380 }
381
382
383 //----------------------------------------------------------------------------
384 void
385 VISU_ColoredPL
386 ::Build() 
387 {
388   myExtractor->SetInput( Superclass::GetClippedInput() );
389   myFieldTransform->SetInput(myExtractor->GetOutput());
390
391   myThreshold->SetInput( myFieldTransform->GetOutput() );
392   // The pass filter is used here for possibility to include/exclude 
393   // threshold filter before it.
394   myPassFilter->SetInput( myFieldTransform->GetOutput() );
395
396   GetMapperHolder()->SetLookupTable(GetMapperTable());
397   //GetMapper()->InterpolateScalarsBeforeMappingOn();
398   GetMapper()->SetUseLookupTableScalarRange( true );
399   GetMapper()->SetColorModeToMapScalars();
400   GetMapper()->ScalarVisibilityOn();
401 }
402
403
404 //----------------------------------------------------------------------------
405 void
406 VISU_ColoredPL
407 ::Update() 
408
409   vtkFloatingPointType *aRange = GetScalarRange();
410   vtkFloatingPointType aScalarRange[2] = {aRange[0], aRange[1]};
411   if(myBarTable->GetScale() == VTK_SCALE_LOG10)
412     VISU_LookupTable::ComputeLogRange(aRange, aScalarRange);
413
414   if(!VISU::CheckIsSameRange(myMapperTable->GetRange(), aScalarRange))
415     myMapperTable->SetRange(aScalarRange);
416   
417   myMapperTable->Build();
418   myBarTable->Build();
419
420   Superclass::Update();
421 }
422
423
424 //----------------------------------------------------------------------------
425 unsigned long int
426 VISU_ColoredPL
427 ::GetMemorySize()
428 {
429   unsigned long int aSize = Superclass::GetMemorySize();
430
431   if(vtkDataObject* aDataObject = myExtractor->GetInput())
432     aSize = aDataObject->GetActualMemorySize() * 1024;
433   
434   if(vtkDataObject* aDataObject = myFieldTransform->GetInput())
435     aSize += aDataObject->GetActualMemorySize() * 1024;
436   
437   return aSize;
438 }
439
440
441 //----------------------------------------------------------------------------
442 VISU_LookupTable *
443 VISU_ColoredPL
444 ::GetMapperTable()
445
446   return myMapperTable.GetPointer();
447 }
448
449
450 //----------------------------------------------------------------------------
451 VISU_LookupTable*
452 VISU_ColoredPL
453 ::GetBarTable()
454 {
455   return myBarTable.GetPointer();
456 }
457
458
459 //----------------------------------------------------------------------------
460 VISU_Extractor*
461 VISU_ColoredPL
462 ::GetExtractorFilter()
463 {
464   return myExtractor.GetPointer();
465 }
466
467
468 //----------------------------------------------------------------------------
469 VISU_FieldTransform*
470 VISU_ColoredPL
471 ::GetFieldTransformFilter()
472 {
473   return myFieldTransform.GetPointer();
474 }
475
476
477 //----------------------------------------------------------------------------
478 void 
479 VISU_ColoredPL
480 ::SetMapScale(vtkFloatingPointType theMapScale)
481 {
482   if(!VISU::CheckIsSameValue(myMapperTable->GetMapScale(), theMapScale)){
483     myMapperTable->SetMapScale(theMapScale);
484     myMapperTable->Build();
485   }
486 }
487
488 vtkFloatingPointType
489 VISU_ColoredPL
490 ::GetMapScale()
491 {
492   return myMapperTable->GetMapScale();
493 }
494
495
496 //----------------------------------------------------------------------------
497 void
498 VISU_ColoredPL
499 ::GetSourceRange(vtkFloatingPointType theRange[2])
500 {
501   myExtractor->Update();
502   myExtractor->GetOutput()->GetScalarRange( theRange );
503   
504   if (isnan(theRange[0]) || isnan(theRange[1]))
505     throw std::runtime_error("Arithmetic exception detected");
506 }
507
508 void
509 VISU_ColoredPL
510 ::SetSourceRange()
511 {
512   vtkFloatingPointType aRange[2];
513   GetSourceRange( aRange );
514   SetScalarRange( aRange );
515 }