Salome HOME
Copyrights update 2015.
[modules/smesh.git] / src / OBJECT / SMESH_ScalarBarActor.cxx
1 // Copyright (C) 2007-2015  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, or (at your option) any later version.
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 //  File   : SMESH_ScalarBarActor.cxx
23 //  Author : Roman NIKOLAEV
24 //  Module : SMESH
25 //
26
27 #include "SMESH_ScalarBarActor.h"
28
29 #include <vtkCellArray.h>
30 #include <vtkCellData.h>
31 #include <vtkObjectFactory.h>
32 #include <vtkPolyData.h>
33 #include <vtkPolyDataMapper2D.h>
34 #include <vtkScalarsToColors.h>
35 #include <vtkTextMapper.h>
36 #include <vtkTextProperty.h>
37 #include <vtkViewport.h>
38 #include <vtkWindow.h>
39 #include <vtkLookupTable.h>
40 #include <vtkProperty2D.h>
41
42 #define SHRINK_COEF 0.08;
43
44 vtkStandardNewMacro(SMESH_ScalarBarActor);
45
46 vtkCxxSetObjectMacro(SMESH_ScalarBarActor,LookupTable,vtkScalarsToColors);
47 vtkCxxSetObjectMacro(SMESH_ScalarBarActor,LabelTextProperty,vtkTextProperty);
48 vtkCxxSetObjectMacro(SMESH_ScalarBarActor,TitleTextProperty,vtkTextProperty);
49
50 //----------------------------------------------------------------------------
51 // Instantiate object with 64 maximum colors; 5 labels; %%-#6.3g label
52 // format, no title, and vertical orientation. The initial scalar bar
53 // size is (0.05 x 0.8) of the viewport size.
54 SMESH_ScalarBarActor::SMESH_ScalarBarActor() {
55   this->LookupTable = NULL;
56   this->Position2Coordinate->SetValue(0.17, 0.8);
57   
58   this->PositionCoordinate->SetCoordinateSystemToNormalizedViewport();
59   this->PositionCoordinate->SetValue(0.82,0.1);
60   
61   this->MaximumNumberOfColors = 64;
62   this->NumberOfLabels = 5;
63   this->NumberOfLabelsBuilt = 0;
64   this->Orientation = VTK_ORIENT_VERTICAL;
65   this->Title = NULL;
66
67   this->LabelTextProperty = vtkTextProperty::New();
68   this->LabelTextProperty->SetFontSize(12);
69   this->LabelTextProperty->SetBold(1);
70   this->LabelTextProperty->SetItalic(1);
71   this->LabelTextProperty->SetShadow(1);
72   this->LabelTextProperty->SetFontFamilyToArial();
73
74   this->TitleTextProperty = vtkTextProperty::New();
75   this->TitleTextProperty->ShallowCopy(this->LabelTextProperty);
76
77   this->LabelFormat = new char[8]; 
78   sprintf(this->LabelFormat,"%s","%-#6.3g");
79
80   this->TitleMapper = vtkTextMapper::New();
81   this->TitleActor = vtkActor2D::New();
82   this->TitleActor->SetMapper(this->TitleMapper);
83   this->TitleActor->GetPositionCoordinate()->
84     SetReferenceCoordinate(this->PositionCoordinate);
85   
86   this->TextMappers = NULL;
87   this->TextActors = NULL;
88
89   this->ScalarBar = vtkPolyData::New();
90   this->ScalarBarMapper = vtkPolyDataMapper2D::New();
91   this->ScalarBarMapper->SetInputData(this->ScalarBar);
92   this->ScalarBarActor = vtkActor2D::New();
93   this->ScalarBarActor->SetMapper(this->ScalarBarMapper);
94   this->ScalarBarActor->GetPositionCoordinate()->
95     SetReferenceCoordinate(this->PositionCoordinate);
96   this->LastOrigin[0] = 0;
97   this->LastOrigin[1] = 0;
98   this->LastSize[0] = 0;
99   this->LastSize[1] = 0;
100
101
102   // rnv begin
103   // Customization of the vtkScalarBarActor to show distribution histogram.
104   myDistribution = vtkPolyData::New();
105   myDistributionMapper = vtkPolyDataMapper2D::New();
106   myDistributionMapper->SetInputData(this->myDistribution);
107   
108   myDistributionActor = vtkActor2D::New();
109   myDistributionActor->SetMapper(this->myDistributionMapper);
110   myDistributionActor->GetPositionCoordinate()->
111     SetReferenceCoordinate(this->PositionCoordinate);
112
113   // By default distribution histogram is invisible
114   myDistributionActor->SetVisibility(0);
115
116   // By default monocolor
117   myDistributionColoringType = SMESH_MONOCOLOR_TYPE;
118   // rnv end
119 }
120
121 //----------------------------------------------------------------------------
122 // Release any graphics resources that are being consumed by this actor.
123 // The parameter window could be used to determine which graphic
124 // resources to release.
125 void SMESH_ScalarBarActor::ReleaseGraphicsResources(vtkWindow *win)
126 {
127   this->TitleActor->ReleaseGraphicsResources(win);
128   if (this->TextMappers != NULL )
129     {
130     for (int i=0; i < this->NumberOfLabelsBuilt; i++)
131       {
132       this->TextActors[i]->ReleaseGraphicsResources(win);
133       }
134     }
135   this->ScalarBarActor->ReleaseGraphicsResources(win);
136   // rnv begin
137   // Customization of the vtkScalarBarActor to show distribution histogram.
138   myDistributionActor->ReleaseGraphicsResources(win);
139 }
140
141
142 /*--------------------------------------------------------------------------*/
143 SMESH_ScalarBarActor::~SMESH_ScalarBarActor() {
144   if (this->LabelFormat) 
145     {
146     delete [] this->LabelFormat;
147     this->LabelFormat = NULL;
148     }
149
150   this->TitleMapper->Delete();
151   this->TitleActor->Delete();
152
153   if (this->TextMappers != NULL )
154     {
155     for (int i=0; i < this->NumberOfLabelsBuilt; i++)
156       {
157       this->TextMappers[i]->Delete();
158       this->TextActors[i]->Delete();
159       }
160     delete [] this->TextMappers;
161     delete [] this->TextActors;
162     }
163
164   this->ScalarBar->Delete();
165   this->ScalarBarMapper->Delete();
166   this->ScalarBarActor->Delete();
167
168   if (this->Title)
169     {
170     delete [] this->Title;
171     this->Title = NULL;
172     }
173   
174   this->SetLookupTable(NULL);
175   this->SetLabelTextProperty(NULL);
176   this->SetTitleTextProperty(NULL);
177   
178   // rnv begin
179   // Customization of the vtkScalarBarActor to show distribution histogram:
180   myDistribution->Delete();
181   myDistributionMapper->Delete();
182   myDistributionActor->Delete();
183   // rnv end
184 }
185
186 //----------------------------------------------------------------------------
187 int SMESH_ScalarBarActor::RenderOverlay(vtkViewport *viewport)
188 {
189   int renderedSomething = 0;
190   int i;
191   
192   // Everything is built, just have to render
193   if (this->Title != NULL)
194     {
195     renderedSomething += this->TitleActor->RenderOverlay(viewport);
196     }
197   this->ScalarBarActor->RenderOverlay(viewport);
198   this->myDistributionActor->RenderOverlay(viewport);
199   if( this->TextActors == NULL)
200     {
201      vtkWarningMacro(<<"Need a mapper to render a scalar bar");
202      return renderedSomething;
203     }
204   
205   for (i=0; i<this->NumberOfLabels; i++)
206     {
207     renderedSomething += this->TextActors[i]->RenderOverlay(viewport);
208     }
209
210   renderedSomething = (renderedSomething > 0)?(1):(0);
211
212   return renderedSomething;
213 }
214
215
216 //----------------------------------------------------------------------------
217 int SMESH_ScalarBarActor::RenderOpaqueGeometry(vtkViewport *viewport)
218 {
219   int renderedSomething = 0;
220   int i;
221   int size[2];
222   
223   if (!this->LookupTable)
224     {
225     vtkWarningMacro(<<"Need a mapper to render a scalar bar");
226     return 0;
227     }
228
229   if (!this->TitleTextProperty)
230     {
231     vtkErrorMacro(<<"Need title text property to render a scalar bar");
232     return 0;
233     }
234
235   if (!this->LabelTextProperty)
236     {
237     vtkErrorMacro(<<"Need label text property to render a scalar bar");
238     return 0;
239     }
240
241   // Check to see whether we have to rebuild everything
242   int positionsHaveChanged = 0;
243   if (viewport->GetMTime() > this->BuildTime || 
244       (viewport->GetVTKWindow() && 
245        viewport->GetVTKWindow()->GetMTime() > this->BuildTime))
246     {
247     // if the viewport has changed we may - or may not need
248     // to rebuild, it depends on if the projected coords chage
249     int *barOrigin;
250     barOrigin = this->PositionCoordinate->GetComputedViewportValue(viewport);
251     size[0] = 
252       this->Position2Coordinate->GetComputedViewportValue(viewport)[0] -
253       barOrigin[0];
254     size[1] = 
255       this->Position2Coordinate->GetComputedViewportValue(viewport)[1] -
256       barOrigin[1];
257     if (this->LastSize[0] != size[0] || 
258         this->LastSize[1] != size[1] ||
259         this->LastOrigin[0] != barOrigin[0] || 
260         this->LastOrigin[1] != barOrigin[1])
261       {
262       positionsHaveChanged = 1;
263       }
264     }
265   
266   // Check to see whether we have to rebuild everything
267   if (positionsHaveChanged ||
268       this->GetMTime() > this->BuildTime || 
269       this->LookupTable->GetMTime() > this->BuildTime ||
270       this->LabelTextProperty->GetMTime() > this->BuildTime ||
271       this->TitleTextProperty->GetMTime() > this->BuildTime)
272     {
273     vtkDebugMacro(<<"Rebuilding subobjects");
274
275     // Delete previously constructed objects
276     //
277     if (this->TextMappers != NULL )
278       {
279       for (i=0; i < this->NumberOfLabelsBuilt; i++)
280         {
281         this->TextMappers[i]->Delete();
282         this->TextActors[i]->Delete();
283         }
284       delete [] this->TextMappers;
285       delete [] this->TextActors;
286       }
287
288     // Build scalar bar object; determine its type
289     //
290     // is this a vtkLookupTable or a subclass of vtkLookupTable 
291     // with its scale set to log
292     // NOTE: it's possible we could to without the 'lut' variable
293     // later in the code, but if the vtkLookupTableSafeDownCast operation
294     // fails for some reason, this code will break in new ways. So, the 'LUT'
295     // variable is used for this operation only
296     vtkLookupTable *LUT = vtkLookupTable::SafeDownCast( this->LookupTable );
297     int isLogTable = 0;
298     if ( LUT )
299       {
300       if ( LUT->GetScale() == VTK_SCALE_LOG10 )
301         {
302         isLogTable = 1; 
303         }
304       }
305     
306     // we hard code how many steps to display
307     vtkScalarsToColors *lut = this->LookupTable;
308     int numColors = this->MaximumNumberOfColors;
309     double *range = lut->GetRange();
310
311     int numPts = 2*(numColors + 1);
312     vtkPoints *pts = vtkPoints::New();
313     pts->SetNumberOfPoints(numPts);
314     vtkCellArray *polys = vtkCellArray::New();
315     polys->Allocate(polys->EstimateSize(numColors,4));
316     vtkUnsignedCharArray *colors = vtkUnsignedCharArray::New();
317     colors->SetNumberOfComponents(3);
318     colors->SetNumberOfTuples(numColors);
319
320
321     // rnv begin
322     // Customization of the vtkScalarBarActor to show distribution histogram.
323     bool distrVisibility =  (numColors == this->myNbValues.size());
324     vtkPoints *distrPts;
325     vtkCellArray *distrPolys;
326     vtkUnsignedCharArray *distColors = 0;
327     int numDistrPts = 0, numPositiveVal=0, maxValue=0;
328     if(!distrVisibility)
329       vtkDebugMacro(<<" Distribution invisible, because numColors == this->myNbValues.size()");
330
331     if (distrVisibility && GetDistributionVisibility()) {
332       for( i=0 ;i<myNbValues.size();i++ ) {
333         if(myNbValues[i]) {
334           numPositiveVal++;
335           maxValue = std::max(maxValue,myNbValues[i]);
336         } 
337       }
338       numDistrPts = 4*(numPositiveVal);
339       distrPts = vtkPoints::New();
340       distrPolys = vtkCellArray::New();
341       distrPts->SetNumberOfPoints(numDistrPts);
342       distrPolys->Allocate(distrPolys->EstimateSize(numPositiveVal,4));
343       this->myDistribution->Initialize();
344       this->myDistribution->SetPoints(distrPts);
345       this->myDistribution->SetPolys(distrPolys);
346       distrPts->Delete();
347       distrPolys->Delete();
348       if ( myDistributionColoringType == SMESH_MULTICOLOR_TYPE ) {
349         distColors = vtkUnsignedCharArray::New();
350         distColors->SetNumberOfComponents(3);
351         distColors->SetNumberOfTuples(numPositiveVal);
352         this->myDistribution->GetCellData()->SetScalars(distColors);
353         distColors->Delete();
354       } else if( myDistributionColoringType == SMESH_MONOCOLOR_TYPE ){
355         this->myDistribution->GetCellData()->SetScalars(NULL);
356       }
357     } else {
358       myDistribution->Reset();
359     }
360     // rnv end
361
362     this->ScalarBarActor->SetProperty(this->GetProperty());
363     this->ScalarBar->Initialize();
364     this->ScalarBar->SetPoints(pts);
365     this->ScalarBar->SetPolys(polys);
366     this->ScalarBar->GetCellData()->SetScalars(colors);
367     pts->Delete(); polys->Delete(); colors->Delete();
368
369     // get the viewport size in display coordinates
370     int *barOrigin, barWidth, barHeight, distrHeight;
371     barOrigin = this->PositionCoordinate->GetComputedViewportValue(viewport);
372     size[0] = 
373       this->Position2Coordinate->GetComputedViewportValue(viewport)[0] -
374       barOrigin[0];
375     size[1] = 
376       this->Position2Coordinate->GetComputedViewportValue(viewport)[1] -
377       barOrigin[1];
378     this->LastOrigin[0] = barOrigin[0];
379     this->LastOrigin[1] = barOrigin[1];
380     this->LastSize[0] = size[0];
381     this->LastSize[1] = size[1];
382     
383     // Update all the composing objects
384     this->TitleActor->SetProperty(this->GetProperty());
385     this->TitleMapper->SetInput(this->Title);
386     if (this->TitleTextProperty->GetMTime() > this->BuildTime)
387       {
388       // Shallow copy here so that the size of the title prop is not affected
389       // by the automatic adjustment of its text mapper's size (i.e. its
390       // mapper's text property is identical except for the font size
391       // which will be modified later). This allows text actors to
392       // share the same text property, and in that case specifically allows
393       // the title and label text prop to be the same.
394       this->TitleMapper->GetTextProperty()->ShallowCopy(this->TitleTextProperty);
395       this->TitleMapper->GetTextProperty()->SetJustificationToCentered();
396       }
397     
398     // find the best size for the title font
399     int titleSize[2];
400     this->SizeTitle(titleSize, size, viewport);
401     
402     // find the best size for the ticks
403     int labelSize[2];
404     this->AllocateAndSizeLabels(labelSize, size, viewport,range);
405     this->NumberOfLabelsBuilt = this->NumberOfLabels;
406     
407     // generate points
408     double x[3]; x[2] = 0.0;
409     double delta, itemH, shrink;
410     if ( this->Orientation == VTK_ORIENT_VERTICAL ) {
411       // rnv begin
412       // Customization of the vtkScalarBarActor to show distribution histogram.
413       double delimeter=0.0;
414       if(GetDistributionVisibility() && distrVisibility) {
415         delimeter=0.01*size[0]; //1 % from horizontal size of the full presentation size.
416         barWidth = size[0] - 4 - labelSize[0];
417         distrHeight = barWidth/2;
418       } else {
419         barWidth = size[0] - 4 - labelSize[0];
420         distrHeight = 0;
421       }
422
423       barHeight = (int)(0.86*size[1]);
424       delta=(double)barHeight/numColors;
425       
426       for ( i=0; i<numPts/2; i++ ) {
427         x[0] = distrHeight+delimeter/2.0;
428         x[1] = i*delta;
429         pts->SetPoint(2*i,x);
430         x[0] = barWidth;
431         pts->SetPoint(2*i+1,x);
432       }
433
434       if(GetDistributionVisibility() && distrVisibility) {
435         // Distribution points 
436         shrink = delta*SHRINK_COEF;
437         vtkIdType distPtsId=0;
438         vtkIdType distPtsIds[4];
439         for(i=0; i<numColors; i++) {
440           if(myNbValues[i]) {
441             itemH = distrHeight*((double)myNbValues[i]/maxValue);
442             
443             if(distrHeight == itemH) 
444               itemH = itemH - delimeter/2;
445
446             x[1] = i*delta+shrink;
447
448             // first point of polygon (quadrangle)
449             x[0] = 0; 
450             distPtsIds[0] = distPtsId;
451             distrPts->SetPoint(distPtsId++,x);
452
453             // second point of polygon (quadrangle)
454             x[0] = itemH;
455             distPtsIds[1] = distPtsId;
456             distrPts->SetPoint(distPtsId++,x);
457
458             x[1] = i*delta+delta-shrink;
459
460             // third point of polygon (quadrangle)
461             x[0] = 0; 
462             distPtsIds[3] = distPtsId;
463             distrPts->SetPoint(distPtsId++,x);
464
465             // fourth point of polygon (quadrangle)
466             x[0] = itemH;
467             distPtsIds[2] = distPtsId;
468             distrPts->SetPoint(distPtsId++,x);
469
470             //Inser Quadrangle
471             distrPolys->InsertNextCell(4,distPtsIds);
472           }
473         }
474       }    
475     }
476     // rnv end
477     else {
478       barWidth = size[0];
479       
480       // rnv begin
481       // Customization of the vtkScalarBarActor to show distribution histogram.
482       double coef1, delimeter=0.0;
483       if(GetDistributionVisibility() && distrVisibility) {
484         coef1=0.62;
485         distrHeight = (int)((coef1/2)*size[1]);
486         //delimeter between distribution diagram and scalar bar 
487         delimeter=0.02*size[1];
488       }
489       else {
490         coef1=0.4;
491         barHeight = (int)(coef1*size[1]);
492         distrHeight = 0;
493       }
494       
495       barHeight = (int)(coef1*size[1]);
496       
497       delta=(double)barWidth/numColors;
498       for (i=0; i<numPts/2; i++) {
499         x[0] = i*delta;
500         x[1] = barHeight;
501         pts->SetPoint(2*i,x);                        
502         x[1] = distrHeight + delimeter;
503         pts->SetPoint(2*i+1,x);
504       }
505       
506       if(GetDistributionVisibility() && distrVisibility) {
507         // Distribution points 
508         shrink = delta*SHRINK_COEF;
509         vtkIdType distPtsId=0;
510         vtkIdType distPtsIds[4];
511         for(i=0; i<numColors; i++) {
512           if(myNbValues[i]) {
513             itemH = distrHeight*((double)myNbValues[i]/maxValue);
514             
515             // first point of polygon (quadrangle)
516             x[0] = i*delta+shrink; 
517             x[1] = 0;
518             distPtsIds[0] = distPtsId;
519             distrPts->SetPoint(distPtsId++,x);
520             
521             // second point of polygon (quadrangle)
522             x[0] = i*delta+shrink; 
523             x[1] = itemH;
524             distPtsIds[3] = distPtsId;
525             distrPts->SetPoint(distPtsId++,x);
526             
527             // third point of polygon (quadrangle)
528             x[0] = i*delta+delta-shrink; 
529             x[1] = 0;
530             distPtsIds[1] = distPtsId;
531             distrPts->SetPoint(distPtsId++,x);
532             
533             // fourth point of polygon (quadrangle)
534             x[0] = i*delta+delta-shrink; 
535             x[1] = itemH;
536             distPtsIds[2] = distPtsId;
537             distrPts->SetPoint(distPtsId++,x);
538             
539             // Add polygon into poly data
540             distrPolys->InsertNextCell(4,distPtsIds);
541           }
542         } 
543       }
544       // rnv end
545     }
546     
547     //polygons & cell colors
548     unsigned char *rgba, *rgb;
549     vtkIdType ptIds[4], dcCount=0;
550     for (i=0; i<numColors; i++)
551       {
552       ptIds[0] = 2*i;
553       ptIds[1] = ptIds[0] + 1;
554       ptIds[2] = ptIds[1] + 2;
555       ptIds[3] = ptIds[0] + 2;
556       polys->InsertNextCell(4,ptIds);
557
558       if ( isLogTable )
559         {
560         double rgbval = log10(range[0]) + 
561           i*(log10(range[1])-log10(range[0]))/(numColors -1);
562         rgba = lut->MapValue(pow(10.0,rgbval));
563         }
564       else
565         {
566         rgba = lut->MapValue(range[0] + (range[1] - range[0])*
567                              ((double)i /(numColors-1.0)));
568         }
569
570       rgb = colors->GetPointer(3*i); //write into array directly
571       rgb[0] = rgba[0];
572       rgb[1] = rgba[1];
573       rgb[2] = rgba[2];
574       
575       // rnv begin
576       // Customization of the vtkScalarBarActor to show distribution histogram.
577       if(myNbValues[i] && myDistributionColoringType == SMESH_MULTICOLOR_TYPE && GetDistributionVisibility() && distrVisibility)
578         {
579           rgb = distColors->GetPointer(3*dcCount); //write into array directly
580           rgb[0] = rgba[0];
581           rgb[1] = rgba[1];
582           rgb[2] = rgba[2];
583           dcCount++;
584         }
585       }
586
587     // Now position everything properly
588     //
589     double val;
590     if (this->Orientation == VTK_ORIENT_VERTICAL)
591       {
592       int sizeTextData[2];
593       
594       // center the title
595       this->TitleActor->SetPosition(size[0]/2, 0.9*size[1]);
596       
597       for (i=0; i < this->NumberOfLabels; i++)
598         {
599         if (this->NumberOfLabels > 1)
600           {
601           val = (double)i/(this->NumberOfLabels-1) *barHeight;
602           }
603         else 
604           {
605           val = 0.5*barHeight;
606           }
607         this->TextMappers[i]->GetSize(viewport,sizeTextData);
608         this->TextMappers[i]->GetTextProperty()->SetJustificationToLeft();
609         this->TextActors[i]->SetPosition(barWidth+3,
610                                          val - sizeTextData[1]/2);
611         }
612       }
613     else
614       {
615       this->TitleActor->SetPosition(size[0]/2, 
616                                     barHeight + labelSize[1] + 0.1*size[1]);
617       for (i=0; i < this->NumberOfLabels; i++)
618         {
619         this->TextMappers[i]->GetTextProperty()->SetJustificationToCentered();
620         if (this->NumberOfLabels > 1)
621           {
622           val = (double)i/(this->NumberOfLabels-1) * barWidth;
623           }
624         else
625           {
626           val = 0.5*barWidth;
627           }
628         this->TextActors[i]->SetPosition(val, barHeight + 0.05*size[1]);
629         }
630       }
631
632     this->BuildTime.Modified();
633     }
634
635   // Everything is built, just have to render
636   if (this->Title != NULL)
637     {
638     renderedSomething += this->TitleActor->RenderOpaqueGeometry(viewport);
639     }
640   this->ScalarBarActor->RenderOpaqueGeometry(viewport);
641   this->myDistributionActor->RenderOpaqueGeometry(viewport);
642   for (i=0; i<this->NumberOfLabels; i++)
643     {
644     renderedSomething += this->TextActors[i]->RenderOpaqueGeometry(viewport);
645     }
646
647   renderedSomething = (renderedSomething > 0)?(1):(0);
648
649   return renderedSomething;
650 }
651
652 //----------------------------------------------------------------------------
653 void SMESH_ScalarBarActor::PrintSelf(ostream& os, vtkIndent indent)
654 {
655   this->Superclass::PrintSelf(os,indent);
656
657   if ( this->LookupTable )
658     {
659     os << indent << "Lookup Table:\n";
660     this->LookupTable->PrintSelf(os,indent.GetNextIndent());
661     }
662   else
663     {
664     os << indent << "Lookup Table: (none)\n";
665     }
666
667   if (this->TitleTextProperty)
668     {
669     os << indent << "Title Text Property:\n";
670     this->TitleTextProperty->PrintSelf(os,indent.GetNextIndent());
671     }
672   else
673     {
674     os << indent << "Title Text Property: (none)\n";
675     }
676
677   if (this->LabelTextProperty)
678     {
679     os << indent << "Label Text Property:\n";
680     this->LabelTextProperty->PrintSelf(os,indent.GetNextIndent());
681     }
682   else
683     {
684     os << indent << "Label Text Property: (none)\n";
685     }
686
687   os << indent << "Title: " << (this->Title ? this->Title : "(none)") << "\n";
688   os << indent << "Maximum Number Of Colors: " 
689      << this->MaximumNumberOfColors << "\n";
690   os << indent << "Number Of Labels: " << this->NumberOfLabels << "\n";
691   os << indent << "Number Of Labels Built: " << this->NumberOfLabelsBuilt << "\n";
692
693   os << indent << "Orientation: ";
694   if ( this->Orientation == VTK_ORIENT_HORIZONTAL )
695     {
696     os << "Horizontal\n";
697     }
698   else
699     {
700     os << "Vertical\n";
701     }
702
703   os << indent << "Label Format: " << this->LabelFormat << "\n";
704 }
705
706 //----------------------------------------------------------------------------
707 void SMESH_ScalarBarActor::ShallowCopy(vtkProp *prop)
708 {
709   SMESH_ScalarBarActor *a = SMESH_ScalarBarActor::SafeDownCast(prop);
710   if ( a != NULL )
711     {
712     this->SetPosition2(a->GetPosition2());
713     this->SetLookupTable(a->GetLookupTable());
714     this->SetMaximumNumberOfColors(a->GetMaximumNumberOfColors());
715     this->SetOrientation(a->GetOrientation());
716     this->SetLabelTextProperty(a->GetLabelTextProperty());
717     this->SetTitleTextProperty(a->GetTitleTextProperty());
718     this->SetLabelFormat(a->GetLabelFormat());
719     this->SetTitle(a->GetTitle());
720     this->GetPositionCoordinate()->SetCoordinateSystem(
721       a->GetPositionCoordinate()->GetCoordinateSystem());    
722     this->GetPositionCoordinate()->SetValue(
723       a->GetPositionCoordinate()->GetValue());
724     this->GetPosition2Coordinate()->SetCoordinateSystem(
725       a->GetPosition2Coordinate()->GetCoordinateSystem());    
726     this->GetPosition2Coordinate()->SetValue(
727       a->GetPosition2Coordinate()->GetValue());
728     }
729
730   // Now do superclass
731   this->vtkActor2D::ShallowCopy(prop);
732 }
733
734 //----------------------------------------------------------------------------
735 void SMESH_ScalarBarActor::AllocateAndSizeLabels(int *labelSize, 
736                                               int *size,
737                                               vtkViewport *viewport,
738                                               double *range)
739 {
740   labelSize[0] = labelSize[1] = 0;
741
742   this->TextMappers = new vtkTextMapper * [this->NumberOfLabels];
743   this->TextActors = new vtkActor2D * [this->NumberOfLabels];
744
745   char string[512];
746
747   double val;
748   int i;
749   
750   // TODO: this should be optimized, maybe by keeping a list of
751   // allocated mappers, in order to avoid creation/destruction of
752   // their underlying text properties (i.e. each time a mapper is
753   // created, text properties are created and shallow-assigned a font size
754   // which value might be "far" from the target font size).
755
756   // is this a vtkLookupTable or a subclass of vtkLookupTable 
757   // with its scale set to log
758   vtkLookupTable *LUT = vtkLookupTable::SafeDownCast( this->LookupTable );
759   int isLogTable = 0;
760   if ( LUT )
761     {
762     if ( LUT->GetScale() == VTK_SCALE_LOG10 )
763       {
764       isLogTable = 1; 
765       }
766     }
767
768   for (i=0; i < this->NumberOfLabels; i++)
769     {
770     this->TextMappers[i] = vtkTextMapper::New();
771
772     if ( isLogTable )
773       {
774       double lval;
775       if (this->NumberOfLabels > 1)
776         {
777         lval = log10(range[0]) + (double)i/(this->NumberOfLabels-1) *
778           (log10(range[1])-log10(range[0]));
779         }
780       else
781         {
782         lval = log10(range[0]) + 0.5*(log10(range[1])-log10(range[0]));
783         }
784       val = pow(10.0,lval);
785       }
786     else
787       {
788       if (this->NumberOfLabels > 1)
789         {
790         val = range[0] + 
791           (double)i/(this->NumberOfLabels-1) * (range[1]-range[0]);
792         }
793       else
794         {
795         val = range[0] + 0.5*(range[1]-range[0]);
796         }
797       }
798
799     sprintf(string, this->LabelFormat, val);
800     this->TextMappers[i]->SetInput(string);
801
802     // Shallow copy here so that the size of the label prop is not affected
803     // by the automatic adjustment of its text mapper's size (i.e. its
804     // mapper's text property is identical except for the font size
805     // which will be modified later). This allows text actors to
806     // share the same text property, and in that case specifically allows
807     // the title and label text prop to be the same.
808     this->TextMappers[i]->GetTextProperty()->ShallowCopy(
809       this->LabelTextProperty);
810
811     this->TextActors[i] = vtkActor2D::New();
812     this->TextActors[i]->SetMapper(this->TextMappers[i]);
813     this->TextActors[i]->SetProperty(this->GetProperty());
814     this->TextActors[i]->GetPositionCoordinate()->
815       SetReferenceCoordinate(this->PositionCoordinate);
816     }
817
818   if (this->NumberOfLabels)
819     {
820     int targetWidth, targetHeight;
821     // rnv begin
822     // Customization of the vtkScalarBarActor to show distribution histogram.
823     bool distrVisibility = this->MaximumNumberOfColors == this->myNbValues.size();
824     double coef;
825     if( GetDistributionVisibility() && distrVisibility )
826       if(this->Orientation == VTK_ORIENT_VERTICAL)
827         coef = 0.4;
828       else 
829         coef = 0.18;
830     else 
831       if(this->Orientation == VTK_ORIENT_VERTICAL)
832         coef = 0.6;
833       else 
834         coef=0.25;
835
836
837     if ( this->Orientation == VTK_ORIENT_VERTICAL )
838       {
839       targetWidth = (int)(coef*size[0]);
840       targetHeight = (int)(0.86*size[1]/this->NumberOfLabels);
841       }
842     else
843       {
844       targetWidth = (int)(size[0]*0.8/this->NumberOfLabels);
845       targetHeight = (int)(coef*size[1]);
846       }
847     // rnv end
848     
849     vtkTextMapper::SetMultipleConstrainedFontSize(viewport, 
850                                                   targetWidth, 
851                                                   targetHeight,
852                                                   this->TextMappers,
853                                                   this->NumberOfLabels,
854                                                   labelSize);
855     }
856 }
857
858 //----------------------------------------------------------------------------
859 void SMESH_ScalarBarActor::SizeTitle(int *titleSize, 
860                                   int *size, 
861                                   vtkViewport *viewport)
862 {
863   titleSize[0] = titleSize[1] = 0;
864
865   if (this->Title == NULL || !strlen(this->Title))
866     {
867     return;
868     }
869
870   int targetWidth, targetHeight;
871   
872   targetWidth = size[0];
873   // rnv begin
874   // Customization of the vtkScalarBarActor to show distribution histogram.
875   bool distrVisibility =  this->MaximumNumberOfColors == this->myNbValues.size();
876   double coef;
877   if( GetDistributionVisibility() && distrVisibility ) 
878     coef=0.18;
879   else 
880     coef=0.25;
881
882   if ( this->Orientation == VTK_ORIENT_VERTICAL )
883     {
884       targetHeight = (int)(0.1*size[1]);
885     }
886   else
887     {
888       targetHeight = (int)(coef*size[1]);
889     }
890
891   this->TitleMapper->SetConstrainedFontSize(
892     viewport, targetWidth, targetHeight);
893
894   this->TitleMapper->GetSize(viewport, titleSize);
895 }
896
897
898 /*--------------------------------------------------------------------------*/
899 void SMESH_ScalarBarActor::SetDistributionVisibility(int flag) {
900   myDistributionActor->SetVisibility(flag);
901   Modified();
902 }
903
904
905 /*--------------------------------------------------------------------------*/
906 int SMESH_ScalarBarActor::GetDistributionVisibility() {
907   return myDistributionActor->GetVisibility();
908 }
909
910
911 void SMESH_ScalarBarActor::SetDistribution(std::vector<int> theNbValues) {
912   myNbValues = theNbValues;
913
914
915
916 void SMESH_ScalarBarActor::SetDistributionColor (double rgb[3]) {
917   myDistributionActor->GetProperty()->SetColor(rgb);
918   Modified();
919 }
920
921 void SMESH_ScalarBarActor::GetDistributionColor (double rgb[3]) {
922   myDistributionActor->GetProperty()->GetColor(rgb);
923 }