Commit bdb8bbf7 authored by Joachim Pouderoux's avatar Joachim Pouderoux
Browse files

Introduce FunctionBagPlot and related filters.

Add a new plot type to draw functional bag plots.

Add vtkExtractFunctionalBagPlot filter to generate the two 2
components columns needed to fee a functional bag plot.

Add a TransposeTable filter to transpose tables, options allow
to create a new column containing initial row names, and use such
a column to name the new columns.

Change-Id: If7fedd9b99038ab7924f5e9b170ed72bd02f9189
parent 5b4532fe
......@@ -25,6 +25,7 @@ set(Module_SRCS
vtkPlot3D.cxx
vtkPlotBag.cxx
vtkPlotBar.cxx
vtkPlotFunctionalBag.cxx
vtkPlotGrid.cxx
vtkPlotHistogram2D.cxx
vtkPlotLine.cxx
......
......@@ -26,6 +26,7 @@ vtk_add_test_cxx(
TestContextUnicode.cxx
TestControlPointsHandleItem.cxx,-E30
TestDiagram.cxx
TestFunctionalBagPlot.cxx
TestHistogram2D.cxx
TestInteractiveChartXYZ.cxx
TestLegendHiddenPlots.cxx
......
/*=========================================================================
Program: Visualization Toolkit
Module: TestFunctionalBagPlot.cxx
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
#include "vtkChartXY.h"
#include "vtkChartLegend.h"
#include "vtkContextScene.h"
#include "vtkContextView.h"
#include "vtkDoubleArray.h"
#include "vtkLookupTable.h"
#include "vtkPlotFunctionalBag.h"
#include "vtkMath.h"
#include "vtkNew.h"
#include "vtkPen.h"
#include "vtkRenderWindow.h"
#include "vtkRenderWindowInteractor.h"
#include "vtkStringArray.h"
#include "vtkTable.h"
#include <sstream>
//----------------------------------------------------------------------------
int TestFunctionalBagPlot(int, char * [])
{
// Creates an input table
const int numCols = 7;
const int numVals = 100;
vtkNew<vtkTable> inputTable;
vtkNew<vtkDoubleArray> arr[numCols];
for (int i = 0; i < numCols; i++)
{
std::stringstream ss;
ss << "Y" << i;
arr[i]->SetName(ss.str().c_str());
arr[i]->SetNumberOfValues(numVals);
for (int j = 0; j < numVals; j++)
{
arr[i]->SetValue(j, (i+1) *
fabs(sin((j * 2.f *vtkMath::Pi()) /
static_cast<float>(numVals))) * j + i * 20);
}
inputTable->AddColumn(arr[i].GetPointer());
}
// Create a X-axis column
vtkNew<vtkDoubleArray> xArr;
xArr->SetName("X");
xArr->SetNumberOfValues(numVals);
for (int j = 0; j < numVals; j++)
{
xArr->SetValue(j, j * 2.0);
}
inputTable->AddColumn(xArr.GetPointer());
// Create the bag columns
vtkNew<vtkDoubleArray> q3Arr;
q3Arr->SetName("Q3");
q3Arr->SetNumberOfComponents(2);
q3Arr->SetNumberOfTuples(numVals);
vtkNew<vtkDoubleArray> q2Arr;
q2Arr->SetName("Q2");
q2Arr->SetNumberOfComponents(2);
q2Arr->SetNumberOfTuples(numVals);
for (int i = 0; i < numVals; i++)
{
double v0, v1;
v0 = arr[1]->GetVariantValue(i).ToFloat();
v1 = arr[5]->GetVariantValue(i).ToFloat();
q3Arr->SetTuple2(i, v0, v1);
v0 = arr[2]->GetVariantValue(i).ToFloat();
v1 = arr[4]->GetVariantValue(i).ToFloat();
q2Arr->SetTuple2(i, v0, v1);
}
inputTable->AddColumn(q3Arr.GetPointer());
inputTable->AddColumn(q2Arr.GetPointer());
// Set up a 2D scene and add an XY chart to it
vtkNew<vtkContextView> view;
view->GetRenderWindow()->SetSize(400, 400);
view->GetRenderWindow()->SetMultiSamples(0);
vtkNew<vtkChartXY> chart;
view->GetScene()->AddItem(chart.GetPointer());
chart->SetShowLegend(true);
chart->GetLegend()->SetHorizontalAlignment(vtkChartLegend::LEFT);
chart->GetLegend()->SetVerticalAlignment(vtkChartLegend::TOP);
// Create the functional bag plots
vtkNew<vtkPlotFunctionalBag> q3Plot;
q3Plot->SetColor(0.5, 0, 0);
q3Plot->SetInputData(inputTable.GetPointer(), "X", "Q3");
chart->AddPlot(q3Plot.GetPointer());
vtkNew<vtkPlotFunctionalBag> q2Plot;
q2Plot->SetColor(1., 0, 0);
q2Plot->SetInputData(inputTable.GetPointer(), "X", "Q2");
chart->AddPlot(q2Plot.GetPointer());
vtkNew<vtkLookupTable> lookup;
lookup->SetNumberOfColors(numCols);
lookup->SetRange(0, numCols-1);
lookup->Build();
for (int j = 0; j < numCols; j++)
{
vtkNew<vtkPlotFunctionalBag> plot;
double rgb[3];
lookup->GetColor(j, rgb);
plot->SetColor(rgb[0], rgb[1], rgb[2]);
plot->SetInputData(inputTable.GetPointer(), "X",
inputTable->GetColumn(j)->GetName());
chart->AddPlot(plot.GetPointer());
}
// Render the scene
view->GetInteractor()->Initialize();
view->GetInteractor()->Start();
return EXIT_SUCCESS;
}
......@@ -52,7 +52,8 @@ public:
POINTS,
BAR,
STACKED,
BAG};
BAG,
FUNCTIONALBAG};
// Description:
// Enum of valid chart action types.
......
......@@ -34,6 +34,7 @@
#include "vtkPlotBar.h"
#include "vtkPlotBag.h"
#include "vtkPlotFunctionalBag.h"
#include "vtkPlotStacked.h"
#include "vtkPlotLine.h"
#include "vtkPlotPoints.h"
......@@ -1039,6 +1040,13 @@ vtkPlot * vtkChartXY::AddPlot(int type)
plot = bar;
break;
}
case FUNCTIONALBAG:
{
vtkPlotFunctionalBag *bag = vtkPlotFunctionalBag::New();
bag->GetBrush()->SetColor(color.GetData());
plot = bag;
break;
}
case STACKED:
{
vtkPlotStacked *stacked = vtkPlotStacked::New();
......
/*=========================================================================
Program: Visualization Toolkit
Module: vtkPlotFunctionalBag.cxx
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
#include "vtkPlotFunctionalBag.h"
#include "vtkAbstractArray.h"
#include "vtkAxis.h"
#include "vtkBrush.h"
#include "vtkContext2D.h"
#include "vtkContextMapper2D.h"
#include "vtkDoubleArray.h"
#include "vtkLookupTable.h"
#include "vtkObjectFactory.h"
#include "vtkPen.h"
#include "vtkPlotLine.h"
#include "vtkPoints2D.h"
#include "vtkRect.h"
#include "vtkScalarsToColors.h"
#include "vtkTable.h"
//-----------------------------------------------------------------------------
vtkStandardNewMacro(vtkPlotFunctionalBag);
//-----------------------------------------------------------------------------
vtkPlotFunctionalBag::vtkPlotFunctionalBag()
{
this->LookupTable = 0;
this->TooltipDefaultLabelFormat = "%l (%x, %y)";
this->LogX = false;
this->LogY = false;
}
//-----------------------------------------------------------------------------
vtkPlotFunctionalBag::~vtkPlotFunctionalBag()
{
if (this->LookupTable)
{
this->LookupTable->UnRegister(this);
}
}
//-----------------------------------------------------------------------------
void vtkPlotFunctionalBag::Update()
{
if (!this->Visible)
{
return;
}
// Check if we have an input
vtkTable *table = this->Data->GetInput();
if (!table)
{
vtkDebugMacro(<< "Update event called with no input table set.");
return;
}
else if(this->Data->GetMTime() > this->BuildTime ||
table->GetMTime() > this->BuildTime ||
(this->LookupTable && this->LookupTable->GetMTime() > this->BuildTime) ||
this->MTime > this->BuildTime)
{
vtkDebugMacro(<< "Updating cached values.");
this->UpdateTableCache(table);
}
else if ((this->XAxis && this->XAxis->GetMTime() > this->BuildTime) ||
(this->YAxis && this->YAxis->GetMaximum() > this->BuildTime))
{
if (this->LogX != this->XAxis->GetLogScale() ||
this->LogY != this->YAxis->GetLogScale())
{
this->UpdateTableCache(table);
}
}
}
//-----------------------------------------------------------------------------
bool vtkPlotFunctionalBag::UpdateTableCache(vtkTable *table)
{
if (!this->LookupTable)
{
this->CreateDefaultLookupTable();
this->LookupTable->SetRange(0, table->GetNumberOfColumns());
this->LookupTable->Build();
}
this->BagPoints->Reset();
vtkDataArray *array[2] = { 0, 0 };
if (!this->GetDataArrays(table, array))
{
this->BuildTime.Modified();
return false;
}
if (array[1]->GetNumberOfComponents() == 1)
{
// The input array has one component, manage it as a line
this->Line->SetInputData(table,
array[0] ? array[0]->GetName() : "", array[1]->GetName());
this->Line->SetUseIndexForXSeries(this->UseIndexForXSeries);
this->Line->SetMarkerStyle(vtkPlotPoints::NONE);
double rgb[3];
this->GetColor(rgb);
this->Line->SetColor(rgb[0], rgb[1], rgb[2]);
this->Line->SetWidth(this->GetWidth());
this->Line->SetPen(this->Pen);
this->Line->SetBrush(this->Brush);
this->Line->Update();
}
else if (array[1]->GetNumberOfComponents() == 2)
{
// The input array has 2 components, this must be a bag
// with {miny,maxy} tuples
vtkDoubleArray* darr = vtkDoubleArray::SafeDownCast(array[1]);
this->LogX = this->XAxis->GetLogScaleActive();
this->LogY = this->YAxis->GetLogScaleActive();
bool xAbs = this->XAxis->GetUnscaledMinimum() < 0.;
bool yAbs = this->YAxis->GetUnscaledMinimum() < 0.;
if (darr)
{
vtkIdType nbRows = array[1]->GetNumberOfTuples();
this->BagPoints->SetNumberOfPoints(2 * nbRows);
for (vtkIdType i = 0; i < nbRows; i++)
{
double y[2];
darr->GetTuple(i, y);
double x = (!this->UseIndexForXSeries && array[0]) ?
array[0]->GetVariantValue(i).ToDouble() : static_cast<double>(i);
if (this->LogX)
{
x = xAbs ? log10(fabs(x)) : log10(x);
}
if (this->LogY)
{
y[0] = yAbs ? log10(fabs(y[0])) : log10(y[0]);
y[1] = yAbs ? log10(fabs(y[1])) : log10(y[1]);
}
this->BagPoints->SetPoint(2 * i, x, y[0]);
this->BagPoints->SetPoint(2 * i + 1, x, y[1]);
}
}
}
this->BuildTime.Modified();
return true;
}
//-----------------------------------------------------------------------------
bool vtkPlotFunctionalBag::GetDataArrays(vtkTable *table, vtkDataArray *array[2])
{
if (!table)
{
return false;
}
// Get the x and y arrays (index 0 and 1 respectively)
array[0] = this->UseIndexForXSeries ?
0 : this->Data->GetInputArrayToProcess(0, table);
array[1] = this->Data->GetInputArrayToProcess(1, table);
if (!array[0] && !this->UseIndexForXSeries)
{
vtkErrorMacro(<< "No X column is set (index 0).");
return false;
}
else if (!array[1])
{
vtkErrorMacro(<< "No Y column is set (index 1).");
return false;
}
else if (!this->UseIndexForXSeries &&
array[0]->GetNumberOfTuples() != array[1]->GetNumberOfTuples())
{
vtkErrorMacro("The x and y columns must have the same number of elements. "
<< array[0]->GetNumberOfTuples() << ", "
<< array[1]->GetNumberOfTuples());
return false;
}
return true;
}
//-----------------------------------------------------------------------------
bool vtkPlotFunctionalBag::Paint(vtkContext2D *painter)
{
// This is where everything should be drawn, or dispatched to other methods.
vtkDebugMacro(<< "Paint event called in vtkPlotFunctionalBag.");
if (!this->Visible)
{
return false;
}
if (this->BagPoints->GetNumberOfPoints() > 0)
{
unsigned char pcolor[4];
double pwidth = this->Pen->GetWidth();
this->Pen->SetWidth(0.);
painter->ApplyPen(this->Pen);
this->Pen->GetColor(pcolor);
this->Brush->SetColor(pcolor);
painter->ApplyBrush(this->Brush);
painter->DrawQuadStrip(this->BagPoints.GetPointer());
this->Pen->SetWidth(pwidth);
}
else
{
this->Line->Paint(painter);
}
return true;
}
//-----------------------------------------------------------------------------
bool vtkPlotFunctionalBag::PaintLegend(vtkContext2D *painter,
const vtkRectf& rect, int index)
{
if (this->BagPoints->GetNumberOfPoints() > 0)
{
vtkNew<vtkPen> blackPen;
blackPen->SetWidth(1.0);
blackPen->SetColor(0, 0, 0, 255);
painter->ApplyPen(blackPen.GetPointer());
painter->ApplyBrush(this->Brush);
painter->DrawRect(rect[0], rect[1], rect[2], rect[3]);
}
else
{
this->Line->PaintLegend(painter, rect, index);
}
return true;
}
//-----------------------------------------------------------------------------
vtkIdType vtkPlotFunctionalBag::GetNearestPoint(const vtkVector2f& point,
const vtkVector2f& tol,
vtkVector2f* loc)
{
if (this->BagPoints->GetNumberOfPoints() == 0)
{
return this->Line->GetNearestPoint(point, tol, loc);
}
return -1;
}
//-----------------------------------------------------------------------------
void vtkPlotFunctionalBag::GetBounds(double bounds[4])
{
if (this->BagPoints->GetNumberOfPoints() > 0)
{
this->BagPoints->GetBounds(bounds);
if (this->LogX)
{
bounds[0] = log10(bounds[0]);
bounds[1] = log10(bounds[1]);
}
if (this->LogY)
{
bounds[2] = log10(bounds[2]);
bounds[3] = log10(bounds[3]);
}
}
else
{
this->Line->GetBounds(bounds);
}
vtkDebugMacro(<< "Bounds: " << bounds[0] << "\t" << bounds[1] << "\t"
<< bounds[2] << "\t" << bounds[3]);
}
//-----------------------------------------------------------------------------
void vtkPlotFunctionalBag::GetUnscaledInputBounds(double bounds[4])
{
if (this->BagPoints->GetNumberOfPoints() > 0)
{
this->BagPoints->GetBounds(bounds);
}
else
{
this->Line->GetUnscaledInputBounds(bounds);
}
vtkDebugMacro(<< "Bounds: " << bounds[0] << "\t" << bounds[1] << "\t"
<< bounds[2] << "\t" << bounds[3]);
}
//-----------------------------------------------------------------------------
void vtkPlotFunctionalBag::PrintSelf(ostream &os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
}
//-----------------------------------------------------------------------------
void vtkPlotFunctionalBag::SetLookupTable(vtkScalarsToColors *lut)
{
if ( this->LookupTable != lut )
{
if ( this->LookupTable)
{
this->LookupTable->UnRegister(this);
}
this->LookupTable = lut;
if (lut)
{
lut->Register(this);
}
this->Modified();
}
}
//-----------------------------------------------------------------------------
vtkScalarsToColors *vtkPlotFunctionalBag::GetLookupTable()
{
if ( this->LookupTable == 0 )
{
this->CreateDefaultLookupTable();
}
return this->LookupTable;
}
//-----------------------------------------------------------------------------
void vtkPlotFunctionalBag::CreateDefaultLookupTable()
{
if ( this->LookupTable)
{
this->LookupTable->UnRegister(this);
}
this->LookupTable = vtkLookupTable::New();
// Consistent Register/UnRegisters.
this->LookupTable->Register(this);
this->LookupTable->Delete();
}
/*=========================================================================
Program: Visualization Toolkit
Module: vtkPlotFunctionalBag.h
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
// .NAME vtkPlotFunctionalBag - Class for drawing an XY line plot or bag
// given two columns from a vtkTable.
//
// .SECTION Description
// Depending on the number of components, this class will draw either
// a line plot (for 1 component column) or, for two components columns,
// a filled polygonal band (the bag) going from the first to the second
// component on the Y-axis along the X-axis. The filter
// vtkExtractFunctionalBagPlot is intended to create such "bag" columns.
//
// .SECTION See Also
// vtkExtractFunctionalBagPlot
#ifndef __vtkPlotFunctionalBag_h
#define __vtkPlotFunctionalBag_h
#include "vtkChartsCoreModule.h" // For export macro
#include "vtkPlot.h"
#include "vtkNew.h" // Needed to hold SP ivars
class vtkDataArray;
class vtkPlotFuntionalBagInternal;
class vtkPlotLine;
class vtkPoints2D;
class vtkScalarsToColors;
class VTKCHARTSCORE_EXPORT vtkPlotFunctionalBag : public vtkPlot
{
public:
vtkTypeMacro(vtkPlotFunctionalBag, vtkPlot);
virtual void PrintSelf(ostream &os, vtkIndent indent);
// Description:
// Creates a functional bag plot object.
static vtkPlotFunctionalBag *New();
// Description:
// Perform any updates to the item that may be necessary before rendering.
// The scene should take care of calling this on all items before their
// Paint function is invoked.
virtual void Update();
// Description:
// Paint event for the plot, called whenever the chart needs to be drawn.
virtual bool Paint(vtkContext2D *painter);
// Description:
// Paint legend event for the plot, called whenever the legend needs the
// plot items symbol/mark/line drawn. A rect is supplied with the lower left