Commit 6f0c543f authored by Ken Martin's avatar Ken Martin Committed by Kitware Robot

Merge topic 'better_sprites'

96b994f1 Fix shadowed var
57daf0bc Some more features and some cleanup reorg.
44e33c4a Cleanup comments
091552c1 Add ability to adjust the table size
edc24bde Improve the capabilities of the vtkPointGaussianMapper
Acked-by: Kitware Robot's avatarKitware Robot <kwrobot@kitware.com>
Tested-by: Scott Wittenburg's avatarScott Wittenburg <scott.wittenburg@kitware.com>
Reviewed-by: Utkarsh Ayachit's avatarUtkarsh Ayachit <utkarsh.ayachit@kitware.com>
Merge-request: !489
parents ae56f168 96b994f1
......@@ -14,15 +14,28 @@
#include "vtkPointGaussianMapper.h"
#include "vtkObjectFactory.h"
#include "vtkPiecewiseFunction.h"
//-----------------------------------------------------------------------------
vtkAbstractObjectFactoryNewMacro(vtkPointGaussianMapper)
vtkCxxSetObjectMacro(vtkPointGaussianMapper, ScaleFunction, vtkPiecewiseFunction);
vtkCxxSetObjectMacro(vtkPointGaussianMapper, ScalarOpacityFunction, vtkPiecewiseFunction);
//-----------------------------------------------------------------------------
vtkPointGaussianMapper::vtkPointGaussianMapper()
{
this->ScaleArray = 0;
this->DefaultRadius = 1.0;
this->OpacityArray = 0;
this->SplatShaderCode = 0;
this->ScaleFunction = 0;
this->ScaleTableSize = 1024;
this->ScalarOpacityFunction = 0;
this->OpacityTableSize = 1024;
this->ScaleFactor = 1.0;
this->Emissive = 1;
}
......@@ -30,6 +43,10 @@ vtkPointGaussianMapper::vtkPointGaussianMapper()
vtkPointGaussianMapper::~vtkPointGaussianMapper()
{
this->SetScaleArray(0);
this->SetOpacityArray(0);
this->SetSplatShaderCode(0);
this->SetScalarOpacityFunction(0);
this->SetScaleFunction(0);
}
//-----------------------------------------------------------------------------
......@@ -38,6 +55,10 @@ void vtkPointGaussianMapper::PrintSelf(ostream& os, vtkIndent indent)
this->Superclass::PrintSelf(os, indent);
os << indent << "Scale Array: " << (this->ScaleArray ? this->ScaleArray : "(none)") << "\n";
os << indent << "Default Radius: " << this->DefaultRadius << "\n";
os << indent << "Opacity Array: " << (this->OpacityArray ? this->OpacityArray : "(none)") << "\n";
os << indent << "SplatShaderCode: " << (this->SplatShaderCode ? this->SplatShaderCode : "(none)") << "\n";
os << indent << "ScaleFactor: " << this->ScaleFactor << "\n";
os << indent << "Emissive: " << this->Emissive << "\n";
os << indent << "OpacityTableSize: " << this->OpacityTableSize << "\n";
os << indent << "ScaleTableSize: " << this->ScaleTableSize << "\n";
}
......@@ -23,6 +23,8 @@
#include "vtkRenderingCoreModule.h" // For export macro
#include "vtkPolyDataMapper.h"
class vtkPiecewiseFunction;
class VTKRENDERINGCORE_EXPORT vtkPointGaussianMapper : public vtkPolyDataMapper
{
public:
......@@ -30,18 +32,32 @@ public:
vtkTypeMacro(vtkPointGaussianMapper, vtkPolyDataMapper)
void PrintSelf(ostream& os, vtkIndent indent);
// Description:
// Set/Get the optional scale transfer function. This is only
// used when a ScaleArray is also specified.
void SetScaleFunction(vtkPiecewiseFunction *);
vtkGetObjectMacro(ScaleFunction,vtkPiecewiseFunction);
// Description:
// The size of the table used in computing scale, used when
// converting a vtkPiecewiseFunction to a table
vtkSetMacro(ScaleTableSize, int);
vtkGetMacro(ScaleTableSize, int);
// Description:
// Convenience method to set the array to scale with.
vtkSetStringMacro(ScaleArray);
vtkGetStringMacro(ScaleArray);
// Description:
// Set the default radius of the point gaussians. This is used if the
// array to scale with has not been set or is set to NULL. If there
// is no scale array and the default radius is set to zero then
// the splats wil be rendered as simple points requiring less memory.
vtkSetMacro(DefaultRadius,double);
vtkGetMacro(DefaultRadius,double);
// Set the default scale factor of the point gaussians. This
// defaults to 1.0. All radius computations will be scaled by the factor
// including the ScaleArray. If a vtkPiecewideFunction is used the
// scaling happens prior to the function lookup.
// A scale factor of 0.0 indicates that the splats should be rendered
// as simple points.
vtkSetMacro(ScaleFactor,double);
vtkGetMacro(ScaleFactor,double);
// Description:
// Treat the points/splats as emissive light sources. The default is true.
......@@ -49,12 +65,50 @@ public:
vtkGetMacro(Emissive, int);
vtkBooleanMacro(Emissive, int);
// Description:
// Set/Get the optional opacity transfer function. This is only
// used when an OpacityArray is also specified.
void SetScalarOpacityFunction(vtkPiecewiseFunction *);
vtkGetObjectMacro(ScalarOpacityFunction,vtkPiecewiseFunction);
// Description:
// The size of the table used in computing opacities, used when
// converting a vtkPiecewiseFunction to a table
vtkSetMacro(OpacityTableSize, int);
vtkGetMacro(OpacityTableSize, int);
// Description:
// Method to set the optional opacity array. If specified this
// array will be used to generate the opacity values.
vtkSetStringMacro(OpacityArray);
vtkGetStringMacro(OpacityArray);
// Description:
// Method to override the fragment shader code for the splat. You can
// set this to draw other shapes. For the OPenGL2 backend some of
// the variables you can use and/or modify include,
// opacity - 0.0 to 1.0
// diffuseColor - vec3
// ambientColor - vec3
// offsetVCVSOutput - vec2 offset in view coordinates from the splat center
vtkSetStringMacro(SplatShaderCode);
vtkGetStringMacro(SplatShaderCode);
protected:
vtkPointGaussianMapper();
~vtkPointGaussianMapper();
char *ScaleArray;
double DefaultRadius;
char *OpacityArray;
char *SplatShaderCode;
vtkPiecewiseFunction *ScaleFunction;
int ScaleTableSize;
vtkPiecewiseFunction *ScalarOpacityFunction;
int OpacityTableSize;
double ScaleFactor;
int Emissive;
private:
......
vtk_add_test_cxx(${vtk-module}CxxTests tests
#TestRenderWidget.cxx # Very experimental, fails, does nothing useful yet.
TestPointGaussianMapper.cxx
TestPointGaussianMapperOpacity.cxx
TestVBOPLYMapper.cxx
TestVBOPointsLines.cxx
TestGaussianBlurPass.cxx
......
......@@ -16,7 +16,8 @@
// .SECTION Thanks
// <verbatim>
//
// This file is part of the PointSprites plugin developed and contributed by
// This file is based loosely on the PointSprites plugin developed
// and contributed by
//
// Copyright (c) CSCS - Swiss National Supercomputing Centre
// EDF - Electricite de France
......@@ -25,8 +26,6 @@
// Stephane Ploix (EDF)
//
// </verbatim>
// .SECTION Description
// this program tests the point sprite support by vtkPointSpriteProperty.
#include "vtkActor.h"
......
/*=========================================================================
Program: Visualization Toolkit
Module: TestSprites.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.
=========================================================================*/
// .SECTION Thanks
// <verbatim>
//
// This file is based loosely on the PointSprites plugin developed
// and contributed by
//
// Copyright (c) CSCS - Swiss National Supercomputing Centre
// EDF - Electricite de France
//
// John Biddiscombe, Ugo Varetto (CSCS)
// Stephane Ploix (EDF)
//
// </verbatim>
#include "vtkActor.h"
#include "vtkCamera.h"
#include "vtkColorTransferFunction.h"
#include "vtkDataObject.h"
#include "vtkDataSetAttributes.h"
#include "vtkNew.h"
#include "vtkPiecewiseFunction.h"
#include "vtkPointGaussianMapper.h"
#include "vtkPointSource.h"
#include "vtkProperty.h"
#include "vtkRandomAttributeGenerator.h"
#include "vtkRenderWindow.h"
#include "vtkRenderWindowInteractor.h"
#include "vtkRenderer.h"
#include "vtkSphereSource.h"
#include "vtkTimerLog.h"
#include "vtkTestUtilities.h"
#include "vtkRegressionTestImage.h"
int TestPointGaussianMapperOpacity(int argc, char *argv[])
{
int desiredPoints = 1.0e4;
vtkNew<vtkPointSource> points;
points->SetNumberOfPoints(desiredPoints);
points->SetRadius(pow(desiredPoints,0.33)*10.0);
points->Update();
vtkNew<vtkRandomAttributeGenerator> randomAttr;
randomAttr->SetInputConnection(points->GetOutputPort());
vtkNew<vtkPointGaussianMapper> mapper;
vtkNew<vtkRenderer> renderer;
renderer->SetBackground(0.0, 0.0, 0.0);
vtkNew<vtkRenderWindow> renderWindow;
renderWindow->SetSize(300, 300);
renderWindow->SetMultiSamples(0);
renderWindow->AddRenderer(renderer.Get());
vtkNew<vtkRenderWindowInteractor> iren;
iren->SetRenderWindow(renderWindow.Get());
vtkNew<vtkActor> actor;
actor->SetMapper(mapper.Get());
renderer->AddActor(actor.Get());
randomAttr->SetDataTypeToFloat();
randomAttr->GeneratePointScalarsOn();
randomAttr->GeneratePointVectorsOn();
randomAttr->GeneratePointArrayOn();
randomAttr->Update();
mapper->SetInputConnection(randomAttr->GetOutputPort());
mapper->SetColorModeToMapScalars();
mapper->SetScalarModeToUsePointFieldData();
mapper->SelectColorArray("RandomPointVectors");
mapper->SetInterpolateScalarsBeforeMapping(0);
mapper->SetScaleArray("RandomPointScalars");
mapper->SetOpacityArray("RandomPointArray");
mapper->EmissiveOff();
// show other shader examples
// the fragment that is rendered is that of a triangle
// large enough to encompass a circle of radius 3
mapper->SetSplatShaderCode(
// this first line keeps the default color opacity calcs
// which you can then modify with additional code below
"//VTK::Color::Impl\n"
// example of a circle with black edges
// " float dist = sqrt(dot(offsetVCVSOutput.xy,offsetVCVSOutput.xy));\n"
// " if (dist > 1.1) { discard; }\n"
// " if (dist < 0.5) { discard; }\n"
// // apply a black edge around the circle
// " if (dist > 1.0 || dist < 0.6) { diffuseColor = vec3(0,0,0); ambientColor = vec3(0,0,0); }\n"
// example for a square
" if (abs(offsetVCVSOutput.x) > 1.0 || abs(offsetVCVSOutput.y) > 1.0) { discard; }\n"
" if (abs(offsetVCVSOutput.x) < 0.6 && abs(offsetVCVSOutput.y) < 0.6) { discard; }\n"
);
vtkNew<vtkColorTransferFunction> ctf;
ctf->AddHSVPoint(0.0,0.1,0.7,1.0);
ctf->AddHSVPoint(1.0,0.9,0.7,1.0);
ctf->SetColorSpaceToHSV();
ctf->HSVWrapOff();
mapper->SetLookupTable(ctf.Get());
vtkNew<vtkPiecewiseFunction> otf;
otf->AddPoint(0.0,0.3);
otf->AddPoint(1.0,1.0);
mapper->SetScalarOpacityFunction(otf.Get());
vtkNew<vtkTimerLog> timer;
timer->StartTimer();
renderWindow->Render();
timer->StopTimer();
double firstRender = timer->GetElapsedTime();
cerr << "first render time: " << firstRender << endl;
timer->StartTimer();
int numRenders = 85;
for (int i = 0; i < numRenders; ++i)
{
renderer->GetActiveCamera()->Azimuth(1);
renderer->GetActiveCamera()->Elevation(1);
renderWindow->Render();
}
timer->StopTimer();
double elapsed = timer->GetElapsedTime();
int numPts = mapper->GetInput()->GetPoints()->GetNumberOfPoints();
cerr << "interactive render time: " << elapsed / numRenders << endl;
cerr << "number of points: " << numPts << endl;
cerr << "points per second: " << numPts*(numRenders/elapsed) << endl;
renderer->GetActiveCamera()->SetPosition(0,0,1);
renderer->GetActiveCamera()->SetFocalPoint(0,0,0);
renderer->GetActiveCamera()->SetViewUp(0,1,0);
renderer->ResetCamera();
// renderer->GetActiveCamera()->Print(cerr);
renderer->GetActiveCamera()->Zoom(10.0);
renderWindow->Render();
int retVal = vtkRegressionTestImage( renderWindow.Get() );
if ( retVal == vtkRegressionTester::DO_INTERACTOR)
{
iren->Start();
}
return !retVal;
}
Markdown is supported
0% or
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment