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