Skip to content
Snippets Groups Projects
Commit e50c7ce8 authored by Nicolas Vuaille's avatar Nicolas Vuaille Committed by Kitware Robot
Browse files

Merge topic '2DGrid'


fdab76cf vktLegendScaleActor expose more Axis parameters
e8193228 vtkAxisActor2D consider UseFontSizeFromProperty for labels
ee72bcd1 Refactor member declaration and initialisation

Acked-by: default avatarKitware Robot <kwrobot@kitware.com>
Tested-by: default avatarbuildbot <buildbot@kitware.com>
Reviewed-by: default avatarThomas Galland <thomas.galland@kitware.com>
Reviewed-by: default avatarCharly Bollinger <charly.bollinger@kitware.com>
Merge-request: !10279
parents 96186a96 fdab76cf
No related branches found
No related tags found
No related merge requests found
## vtkAxisActor2D consider UseFontSizeFromProperty for labels
This property of the vtkAxisActor2D class was used for the title display.
It is now also used for the labels.
......@@ -13,6 +13,7 @@
=========================================================================*/
// This tests the terrain annotation capabilities in VTK.
#include "vtkAxisActor2D.h"
#include "vtkCamera.h"
#include "vtkInteractorStyleTrackballCamera.h"
#include "vtkLegendScaleActor.h"
......@@ -23,39 +24,55 @@
#include "vtkRenderer.h"
#include "vtkSphereSource.h"
#include "vtkTestUtilities.h"
#include "vtkTextProperty.h"
//------------------------------------------------------------------------------
int TestLegendScaleActor(int argc, char* argv[])
{
// Create the RenderWindow, Renderer and both Actors
//
vtkRenderer* ren1 = vtkRenderer::New();
vtkRenderWindow* renWin = vtkRenderWindow::New();
vtkNew<vtkRenderer> ren1;
vtkNew<vtkRenderWindow> renWin;
renWin->SetMultiSamples(0);
renWin->AddRenderer(ren1);
ren1->GetActiveCamera()->ParallelProjectionOn();
vtkInteractorStyleTrackballCamera* style = vtkInteractorStyleTrackballCamera::New();
vtkRenderWindowInteractor* iren = vtkRenderWindowInteractor::New();
vtkNew<vtkInteractorStyleTrackballCamera> style;
vtkNew<vtkRenderWindowInteractor> iren;
iren->SetRenderWindow(renWin);
iren->SetInteractorStyle(style);
// Create a test pipeline
//
vtkSphereSource* ss = vtkSphereSource::New();
vtkPolyDataMapper* mapper = vtkPolyDataMapper::New();
vtkNew<vtkSphereSource> ss;
vtkNew<vtkPolyDataMapper> mapper;
mapper->SetInputConnection(ss->GetOutputPort());
vtkActor* sph = vtkActor::New();
vtkNew<vtkActor> sph;
sph->SetMapper(mapper);
// Create the actor
vtkLegendScaleActor* actor = vtkLegendScaleActor::New();
actor->TopAxisVisibilityOn();
vtkNew<vtkLegendScaleActor> legendActor;
legendActor->TopAxisVisibilityOn();
legendActor->SetLabelModeToXYCoordinates();
legendActor->AllAxesOff();
legendActor->LeftAxisVisibilityOn();
legendActor->TopAxisVisibilityOn();
legendActor->SetLeftBorderOffset(70);
legendActor->SetTopBorderOffset(50);
legendActor->GetTopAxis()->SetNumberOfLabels(3);
legendActor->SetCornerOffsetFactor(1);
vtkNew<vtkTextProperty> textProp;
textProp->SetColor(1, 0.5, 0);
textProp->SetFontSize(18);
textProp->BoldOn();
legendActor->SetUseFontSizeFromProperty(true);
legendActor->SetAxesTextProperty(textProp);
// Add the actors to the renderer, set the background and size
ren1->AddActor(sph);
ren1->AddViewProp(actor);
ren1->AddViewProp(legendActor);
ren1->SetBackground(0.1, 0.2, 0.4);
renWin->SetSize(300, 300);
......@@ -70,14 +87,5 @@ int TestLegendScaleActor(int argc, char* argv[])
iren->Start();
}
ss->Delete();
mapper->Delete();
sph->Delete();
actor->Delete();
style->Delete();
iren->Delete();
renWin->Delete();
ren1->Delete();
return !retVal;
}
b996ef7706b6aaa05403059e7e381552cde10086c7f6f937b0f79070445b4f59b0df6dbf9161e7dc0c14e976cd6f44ac83b01f7e0037ae864f4b543c7e244f69
2824415ba24e1d451774e2faa1248235c57b96a5cfbc07432fbffd3db338a18759386fcd7792b43f70113bce9c8840d9f1783861d86e6ae456f2c8cedc76e2ca
......@@ -44,31 +44,7 @@ vtkAxisActor2D::vtkAxisActor2D()
this->Position2Coordinate->SetValue(0.75, 0.0);
this->Position2Coordinate->SetReferenceCoordinate(nullptr);
this->NumberOfLabels = 5;
this->Title = nullptr;
this->TitlePosition = 0.5;
this->AdjustLabels = 1;
this->TickLength = 5;
this->MinorTickLength = 3;
this->TickOffset = 2;
this->NumberOfMinorTicks = 0;
this->Range[0] = 0.0;
this->Range[1] = 1.0;
this->FontFactor = 1.0;
this->LabelFactor = 0.75;
this->SizeFontRelativeToAxis = 0;
this->UseFontSizeFromProperty = 0;
this->RulerMode = 0;
this->RulerDistance = 1.0;
this->LabelTextProperty = vtkTextProperty::New();
this->LabelTextProperty->SetBold(1);
this->LabelTextProperty->SetItalic(1);
......@@ -86,7 +62,6 @@ vtkAxisActor2D::vtkAxisActor2D()
this->TitleActor->SetMapper(this->TitleMapper);
// To avoid deleting/rebuilding create once up front
this->NumberOfLabelsBuilt = 0;
this->LabelMappers = new vtkTextMapper*[VTK_MAX_LABELS];
this->LabelActors = new vtkActor2D*[VTK_MAX_LABELS];
for (int i = 0; i < VTK_MAX_LABELS; i++)
......@@ -96,22 +71,8 @@ vtkAxisActor2D::vtkAxisActor2D()
this->LabelActors[i]->SetMapper(this->LabelMappers[i]);
}
this->Axis = vtkPolyData::New();
this->AxisMapper = vtkPolyDataMapper2D::New();
this->AxisMapper->SetInputData(this->Axis);
this->AxisActor = vtkActor2D::New();
this->AxisActor->SetMapper(this->AxisMapper);
this->AxisVisibility = 1;
this->TickVisibility = 1;
this->LabelVisibility = 1;
this->TitleVisibility = 1;
this->LastPosition[0] = this->LastPosition[1] = 0;
this->LastPosition2[0] = this->LastPosition2[1] = 0;
this->LastSize[0] = this->LastSize[1] = 0;
this->LastMaxLabelSize[0] = this->LastMaxLabelSize[1] = 0;
}
//------------------------------------------------------------------------------
......@@ -137,10 +98,6 @@ vtkAxisActor2D::~vtkAxisActor2D()
delete[] this->LabelActors;
}
this->Axis->Delete();
this->AxisMapper->Delete();
this->AxisActor->Delete();
this->SetLabelTextProperty(nullptr);
this->SetTitleTextProperty(nullptr);
}
......@@ -537,31 +494,38 @@ void vtkAxisActor2D::BuildAxis(vtkViewport* viewport)
if (positionsHaveChanged || viewportSizeHasChanged ||
this->LabelTextProperty->GetMTime() > this->BuildTime || labeltime > this->BuildTime)
{
if (!this->SizeFontRelativeToAxis)
{
vtkTextMapper::SetMultipleRelativeFontSize(viewport, this->LabelMappers,
this->AdjustedNumberOfLabels, size, this->LastMaxLabelSize,
0.015 * this->FontFactor * this->LabelFactor);
}
else
if (!this->UseFontSizeFromProperty)
{
int minFontSize = 1000, fontSize, minLabel = 0;
for (i = 0; i < this->AdjustedNumberOfLabels; i++)
if (!this->SizeFontRelativeToAxis)
{
fontSize = this->LabelMappers[i]->SetConstrainedFontSize(viewport,
static_cast<int>((1.0 / this->AdjustedNumberOfLabels) * len),
static_cast<int>(0.2 * len));
if (fontSize < minFontSize)
{
minFontSize = fontSize;
minLabel = i;
}
vtkTextMapper::SetMultipleRelativeFontSize(viewport, this->LabelMappers,
this->AdjustedNumberOfLabels, size, this->LastMaxLabelSize,
0.015 * this->FontFactor * this->LabelFactor);
}
for (i = 0; i < this->AdjustedNumberOfLabels; i++)
else
{
this->LabelMappers[i]->GetTextProperty()->SetFontSize(minFontSize);
int minFontSize = 1000, fontSize, minLabel = 0;
for (i = 0; i < this->AdjustedNumberOfLabels; i++)
{
fontSize = this->LabelMappers[i]->SetConstrainedFontSize(viewport,
static_cast<int>((1.0 / this->AdjustedNumberOfLabels) * len),
static_cast<int>(0.2 * len));
if (fontSize < minFontSize)
{
minFontSize = fontSize;
minLabel = i;
}
}
for (i = 0; i < this->AdjustedNumberOfLabels; i++)
{
this->LabelMappers[i]->GetTextProperty()->SetFontSize(minFontSize);
}
this->LabelMappers[minLabel]->GetSize(viewport, this->LastMaxLabelSize);
}
this->LabelMappers[minLabel]->GetSize(viewport, this->LastMaxLabelSize);
}
else
{
this->LabelMappers[0]->GetSize(viewport, this->LastMaxLabelSize);
}
}
......
......@@ -61,6 +61,8 @@
#include "vtkActor2D.h"
#include "vtkRenderingAnnotationModule.h" // For export macro
#include "vtkNew.h" // for vtkNew
VTK_ABI_NAMESPACE_BEGIN
class vtkPolyDataMapper2D;
class vtkPolyData;
......@@ -390,37 +392,37 @@ protected:
vtkTextProperty* LabelTextProperty;
char* Title;
double Range[2];
double TitlePosition;
vtkTypeBool RulerMode;
double RulerDistance;
int NumberOfLabels;
char* LabelFormat;
vtkTypeBool AdjustLabels;
double FontFactor;
double LabelFactor;
int TickLength;
int MinorTickLength;
int TickOffset;
int NumberOfMinorTicks;
double Range[2] = { 0., 1. };
double TitlePosition = 0.5;
vtkTypeBool RulerMode = 0;
double RulerDistance = 1.;
int NumberOfLabels = 5;
vtkTypeBool AdjustLabels = 1;
double FontFactor = 1.;
double LabelFactor = 0.75;
int TickLength = 5;
int MinorTickLength = 3;
int TickOffset = 2;
int NumberOfMinorTicks = 0;
double AdjustedRange[2];
int AdjustedNumberOfLabels;
int NumberOfLabelsBuilt;
int AdjustedNumberOfLabels = 5;
int NumberOfLabelsBuilt = 0;
vtkTypeBool AxisVisibility;
vtkTypeBool TickVisibility;
vtkTypeBool LabelVisibility;
vtkTypeBool TitleVisibility;
vtkTypeBool AxisVisibility = 1;
vtkTypeBool TickVisibility = 1;
vtkTypeBool LabelVisibility = 1;
vtkTypeBool TitleVisibility = 1;
int LastPosition[2];
int LastPosition2[2];
int LastPosition[2] = { 0, 0 };
int LastPosition2[2] = { 0, 0 };
int LastSize[2];
int LastMaxLabelSize[2];
int LastSize[2] = { 0, 0 };
int LastMaxLabelSize[2] = { 0, 0 };
int SizeFontRelativeToAxis;
vtkTypeBool UseFontSizeFromProperty;
int SizeFontRelativeToAxis = 0;
vtkTypeBool UseFontSizeFromProperty = 0;
virtual void BuildAxis(vtkViewport* viewport);
static double ComputeStringOffset(double width, double height, double theta);
......@@ -434,9 +436,9 @@ protected:
vtkTextMapper** LabelMappers;
vtkActor2D** LabelActors;
vtkPolyData* Axis;
vtkPolyDataMapper2D* AxisMapper;
vtkActor2D* AxisActor;
vtkNew<vtkPolyData> Axis;
vtkNew<vtkPolyDataMapper2D> AxisMapper;
vtkNew<vtkActor2D> AxisActor;
vtkTimeStamp AdjustedRangeBuildTime;
vtkTimeStamp BuildTime;
......
......@@ -37,15 +37,6 @@ vtkStandardNewMacro(vtkLegendScaleActor);
//------------------------------------------------------------------------------
vtkLegendScaleActor::vtkLegendScaleActor()
{
this->LabelMode = DISTANCE;
this->RightBorderOffset = 50;
this->TopBorderOffset = 30;
this->LeftBorderOffset = 50;
this->BottomBorderOffset = 30;
this->CornerOffsetFactor = 2.0;
this->RightAxis = vtkAxisActor2D::New();
this->RightAxis->GetPositionCoordinate()->SetCoordinateSystemToViewport();
this->RightAxis->GetPosition2Coordinate()->SetCoordinateSystemToViewport();
this->RightAxis->GetPositionCoordinate()->SetReferenceCoordinate(nullptr);
......@@ -53,7 +44,6 @@ vtkLegendScaleActor::vtkLegendScaleActor()
this->RightAxis->SetNumberOfLabels(5);
this->RightAxis->AdjustLabelsOff();
this->TopAxis = vtkAxisActor2D::New();
this->TopAxis->GetPositionCoordinate()->SetCoordinateSystemToViewport();
this->TopAxis->GetPosition2Coordinate()->SetCoordinateSystemToViewport();
this->TopAxis->GetPositionCoordinate()->SetReferenceCoordinate(nullptr);
......@@ -61,7 +51,6 @@ vtkLegendScaleActor::vtkLegendScaleActor()
this->TopAxis->SetNumberOfLabels(5);
this->TopAxis->AdjustLabelsOff();
this->LeftAxis = vtkAxisActor2D::New();
this->LeftAxis->GetPositionCoordinate()->SetCoordinateSystemToViewport();
this->LeftAxis->GetPosition2Coordinate()->SetCoordinateSystemToViewport();
this->LeftAxis->GetPositionCoordinate()->SetReferenceCoordinate(nullptr);
......@@ -69,7 +58,6 @@ vtkLegendScaleActor::vtkLegendScaleActor()
this->LeftAxis->SetNumberOfLabels(5);
this->LeftAxis->AdjustLabelsOff();
this->BottomAxis = vtkAxisActor2D::New();
this->BottomAxis->GetPositionCoordinate()->SetCoordinateSystemToViewport();
this->BottomAxis->GetPosition2Coordinate()->SetCoordinateSystemToViewport();
this->BottomAxis->GetPositionCoordinate()->SetReferenceCoordinate(nullptr);
......@@ -77,24 +65,14 @@ vtkLegendScaleActor::vtkLegendScaleActor()
this->BottomAxis->SetNumberOfLabels(5);
this->BottomAxis->AdjustLabelsOff();
this->RightAxisVisibility = 1;
this->TopAxisVisibility = 1;
this->LeftAxisVisibility = 1;
this->BottomAxisVisibility = 1;
this->LegendVisibility = 1;
this->Legend = vtkPolyData::New();
this->LegendPoints = vtkPoints::New();
this->Legend->SetPoints(this->LegendPoints);
this->LegendMapper = vtkPolyDataMapper2D::New();
this->LegendMapper->SetInputData(this->Legend);
this->LegendActor = vtkActor2D::New();
this->LegendActor->SetMapper(this->LegendMapper);
// Create the legend
vtkIdType pts[4];
this->LegendPoints->SetNumberOfPoints(10);
vtkCellArray* legendPolys = vtkCellArray::New();
vtkNew<vtkCellArray> legendPolys;
legendPolys->AllocateEstimate(4, 4);
pts[0] = 0;
pts[1] = 1;
......@@ -117,7 +95,6 @@ vtkLegendScaleActor::vtkLegendScaleActor()
pts[3] = 8;
legendPolys->InsertNextCell(4, pts);
this->Legend->SetPolys(legendPolys);
legendPolys->Delete();
// Create the cell data
vtkUnsignedCharArray* colors = vtkUnsignedCharArray::New();
......@@ -131,7 +108,6 @@ vtkLegendScaleActor::vtkLegendScaleActor()
colors->Delete();
// Now the text. The first five are for the 0,1/4,1/2,3/4,1 labels.
this->LegendTitleProperty = vtkTextProperty::New();
this->LegendTitleProperty->SetJustificationToCentered();
this->LegendTitleProperty->SetVerticalJustificationToBottom();
this->LegendTitleProperty->SetBold(1);
......@@ -139,7 +115,7 @@ vtkLegendScaleActor::vtkLegendScaleActor()
this->LegendTitleProperty->SetShadow(1);
this->LegendTitleProperty->SetFontFamilyToArial();
this->LegendTitleProperty->SetFontSize(10);
this->LegendLabelProperty = vtkTextProperty::New();
this->LegendLabelProperty->SetJustificationToCentered();
this->LegendLabelProperty->SetVerticalJustificationToTop();
this->LegendLabelProperty->SetBold(1);
......@@ -161,31 +137,51 @@ vtkLegendScaleActor::vtkLegendScaleActor()
this->LabelMappers[3]->SetInput("3/4");
this->LabelMappers[4]->SetInput("1");
this->Coordinate = vtkCoordinate::New();
this->Coordinate->SetCoordinateSystemToDisplay();
}
//------------------------------------------------------------------------------
vtkLegendScaleActor::~vtkLegendScaleActor()
{
this->RightAxis->Delete();
this->TopAxis->Delete();
this->LeftAxis->Delete();
this->BottomAxis->Delete();
this->Legend->Delete();
this->LegendPoints->Delete();
this->LegendMapper->Delete();
this->LegendActor->Delete();
for (int i = 0; i < 6; i++)
{
this->LabelMappers[i]->Delete();
this->LabelActors[i]->Delete();
}
this->LegendTitleProperty->Delete();
this->LegendLabelProperty->Delete();
this->Coordinate->Delete();
}
//------------------------------------------------------------------------------
void vtkLegendScaleActor::SetAdjustLabels(bool adjust)
{
this->RightAxis->SetAdjustLabels(adjust);
this->TopAxis->SetAdjustLabels(adjust);
this->LeftAxis->SetAdjustLabels(adjust);
this->BottomAxis->SetAdjustLabels(adjust);
}
//------------------------------------------------------------------------------
void vtkLegendScaleActor::SetUseFontSizeFromProperty(bool fromProp)
{
this->RightAxis->SetUseFontSizeFromProperty(fromProp);
this->TopAxis->SetUseFontSizeFromProperty(fromProp);
this->LeftAxis->SetUseFontSizeFromProperty(fromProp);
this->BottomAxis->SetUseFontSizeFromProperty(fromProp);
}
//------------------------------------------------------------------------------
void vtkLegendScaleActor::SetAxesTextProperty(vtkTextProperty* prop)
{
this->RightAxis->SetLabelTextProperty(prop);
this->TopAxis->SetLabelTextProperty(prop);
this->LeftAxis->SetLabelTextProperty(prop);
this->BottomAxis->SetLabelTextProperty(prop);
this->RightAxis->SetTitleTextProperty(prop);
this->TopAxis->SetTitleTextProperty(prop);
this->LeftAxis->SetTitleTextProperty(prop);
this->BottomAxis->SetTitleTextProperty(prop);
this->Modified();
}
//------------------------------------------------------------------------------
......
......@@ -41,6 +41,8 @@
#include "vtkProp.h"
#include "vtkRenderingAnnotationModule.h" // For export macro
#include "vtkNew.h" // for vtkNew
VTK_ABI_NAMESPACE_BEGIN
class vtkAxisActor2D;
class vtkTextProperty;
......@@ -188,6 +190,20 @@ public:
vtkGetObjectMacro(LegendLabelProperty, vtkTextProperty);
///@}
/**
* Configuration forwarded to each axis.
*/
///@{
/// Set the axes text properties.
void SetAxesTextProperty(vtkTextProperty* property);
/// Set the axes to get font size from text property.
void SetUseFontSizeFromProperty(bool sizeFromProp);
/// Set the axes to adjust labels position to a "nice" one.
void SetAdjustLabels(bool ajust);
///@}
///@{
/**
* These are methods to retrieve the vtkAxisActors used to represent
......@@ -215,36 +231,36 @@ protected:
vtkLegendScaleActor();
~vtkLegendScaleActor() override;
int LabelMode;
int RightBorderOffset;
int TopBorderOffset;
int LeftBorderOffset;
int BottomBorderOffset;
double CornerOffsetFactor;
int LabelMode = DISTANCE;
int RightBorderOffset = 50;
int TopBorderOffset = 30;
int LeftBorderOffset = 50;
int BottomBorderOffset = 30;
double CornerOffsetFactor = 2.;
// The four axes around the borders of the renderer
vtkAxisActor2D* RightAxis;
vtkAxisActor2D* TopAxis;
vtkAxisActor2D* LeftAxis;
vtkAxisActor2D* BottomAxis;
vtkNew<vtkAxisActor2D> RightAxis;
vtkNew<vtkAxisActor2D> TopAxis;
vtkNew<vtkAxisActor2D> LeftAxis;
vtkNew<vtkAxisActor2D> BottomAxis;
// Control the display of the axes
vtkTypeBool RightAxisVisibility;
vtkTypeBool TopAxisVisibility;
vtkTypeBool LeftAxisVisibility;
vtkTypeBool BottomAxisVisibility;
vtkTypeBool RightAxisVisibility = 1;
vtkTypeBool TopAxisVisibility = 1;
vtkTypeBool LeftAxisVisibility = 1;
vtkTypeBool BottomAxisVisibility = 1;
// Support for the legend.
vtkTypeBool LegendVisibility;
vtkPolyData* Legend;
vtkPoints* LegendPoints;
vtkPolyDataMapper2D* LegendMapper;
vtkActor2D* LegendActor;
vtkTypeBool LegendVisibility = 1;
vtkNew<vtkPolyData> Legend;
vtkNew<vtkPoints> LegendPoints;
vtkNew<vtkPolyDataMapper2D> LegendMapper;
vtkNew<vtkActor2D> LegendActor;
vtkTextMapper* LabelMappers[6];
vtkActor2D* LabelActors[6];
vtkTextProperty* LegendTitleProperty;
vtkTextProperty* LegendLabelProperty;
vtkCoordinate* Coordinate;
vtkNew<vtkTextProperty> LegendTitleProperty;
vtkNew<vtkTextProperty> LegendLabelProperty;
vtkNew<vtkCoordinate> Coordinate;
vtkTimeStamp BuildTime;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment