Updates will be applied April 15th at 12pm EDT (UTC-0400). GitLab could be a little slow between 12 - 12:45pm EDT.

Commit cd98500c authored by Ken Martin's avatar Ken Martin

add a poly point source where you can explicitly set points

Similar to PolyLineSource which was made a subclass of this.
Useful to have a gui driven point source with explcit points.
parent 68f35cda
......@@ -23,6 +23,7 @@ set(classes
vtkPlatonicSolidSource
vtkPointSource
vtkPolyLineSource
vtkPolyPointSource
vtkProgrammableDataObjectSource
vtkProgrammableSource
vtkRandomHyperTreeGridSource
......
......@@ -18,6 +18,7 @@ vtk_add_test_cxx(vtkFiltersSourcesCxxTests tests
TestPlatonicSolidSource.cxx,NO_VALID
TestPointSource.cxx,NO_VALID
TestPolyLineSource.cxx,NO_VALID
TestPolyPointSource.cxx,NO_VALID
TestProgrammableSource.cxx,NO_VALID
TestRandomHyperTreeGridSource.cxx
TestRectangularButtonSource.cxx,NO_VALID
......
/*=========================================================================
Program: Visualization Toolkit
Module: TestPolyPointSource.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 <cstdio>
#include <vtkPolyPointSource.h>
#include <vtkSmartPointer.h>
int TestPolyPointSource(int vtkNotUsed(argc), char* vtkNotUsed(argv)[])
{
vtkSmartPointer<vtkPolyPointSource> source = vtkSmartPointer<vtkPolyPointSource>::New();
// Basic tests
int expectedNumberOfPoints = 4;
source->SetNumberOfPoints( 4 );
int actualNumberOfPoints = source->GetNumberOfPoints();
if (expectedNumberOfPoints != actualNumberOfPoints)
{
std::cerr << "Expected NumberOfPoints setting to be " << expectedNumberOfPoints << ", got "
<< actualNumberOfPoints << std::endl;
return EXIT_FAILURE;
}
// Test setting individual points
double pts[4][3] = {
{ 1.0, 2.0, 3.0 },
{ 4.0, 5.0, 6.0 },
{ 7.0, 8.0, 9.0 },
{ 10.0, 11.0, 12.0 }
};
for ( int i = 0; i < 4; ++i )
{
source->SetPoint( i, pts[i][0], pts[i][1], pts[i][2] );
}
// Test getting the points
vtkPoints* testPoints = source->GetPoints();
if ( testPoints->GetNumberOfPoints() != 4 )
{
std::cerr << "Expected 4 points in vtkPoints returned from GetPoints() method, but got "
<< testPoints->GetNumberOfPoints() << " instead." << std::endl;
return EXIT_FAILURE;
}
// Test the point values
vtkPoints* outputPoints = source->GetPoints();
for ( int i = 0; i < 4; ++i )
{
double pt[3];
outputPoints->GetPoint(i, pt);
if ( pt[0] != pts[i][0] || pt[1] != pts[i][1] || pt[2] != pts[i][2] )
{
std::cerr << "Point disagreement in point " << i << std::endl;
std::cerr << "Expected point: " << pts[i][0] << ", " << pts[i][1] << ", " << pts[i][2] << std::endl;
std::cerr << "Actual point: " << pt[0] << ", " << pt[1] << ", " << pt[2] << std::endl;
return EXIT_FAILURE;
}
}
// Test setting the points from a vtkPoints object
double newPts[3][3] = {
{ 13.0, 14.0, 15.0 },
{ 16.0, 17.0, 18.0 },
{ 19.0, 20.0, 21.0 }
};
vtkSmartPointer<vtkPoints> newPoints = vtkSmartPointer<vtkPoints>::New();
newPoints->SetNumberOfPoints( 3 );
for ( int i = 0; i < 3; ++i )
{
newPoints->SetPoint( i, newPts[i] );
}
source->SetPoints( newPoints );
actualNumberOfPoints = source->GetNumberOfPoints();
if ( source->GetNumberOfPoints() != 3)
{
std::cerr << "Expected 3 points, got " << source->GetNumberOfPoints() << std::endl;
return EXIT_FAILURE;
}
outputPoints = source->GetPoints();
for ( int i = 0; i < 3; ++i )
{
double pt[3];
outputPoints->GetPoint(i, pt);
if ( pt[0] != newPts[i][0] || pt[1] != newPts[i][1] || pt[2] != newPts[i][2] )
{
std::cerr << "Point disagreement in point " << i << std::endl;
std::cerr << "Expected point: " << newPts[i][0] << ", " << newPts[i][1] << ", " << newPts[i][2] << std::endl;
std::cerr << "Actual point: " << pt[0] << ", " << pt[1] << ", " << pt[2] << std::endl;
return EXIT_FAILURE;
}
}
return EXIT_SUCCESS;
}
......@@ -27,99 +27,12 @@ vtkStandardNewMacro(vtkPolyLineSource);
//----------------------------------------------------------------------------
vtkPolyLineSource::vtkPolyLineSource()
{
this->Points = nullptr;
this->Closed = 0;
this->SetNumberOfInputPorts(0);
}
//----------------------------------------------------------------------------
vtkPolyLineSource::~vtkPolyLineSource()
{
if (this->Points)
{
this->Points->Delete();
}
}
//----------------------------------------------------------------------------
void vtkPolyLineSource::SetNumberOfPoints(vtkIdType numPoints)
{
if (!this->Points)
{
vtkPoints* pts = vtkPoints::New(VTK_DOUBLE);
this->SetPoints(pts);
this->Points = pts;
pts->Delete();
}
if (numPoints != this->GetNumberOfPoints())
{
this->Points->SetNumberOfPoints(numPoints);
this->Modified();
}
}
//----------------------------------------------------------------------------
vtkIdType vtkPolyLineSource::GetNumberOfPoints()
{
if (this->Points)
{
return this->Points->GetNumberOfPoints();
}
return 0;
}
//----------------------------------------------------------------------------
void vtkPolyLineSource::Resize(vtkIdType numPoints)
{
if (!this->Points)
{
this->SetNumberOfPoints(numPoints);
}
if (numPoints != this->GetNumberOfPoints())
{
this->Points->Resize(numPoints);
this->Modified();
}
}
//----------------------------------------------------------------------------
void vtkPolyLineSource::SetPoint(vtkIdType id, double x, double y, double z)
{
if (!this->Points)
{
return;
}
if (id >= this->Points->GetNumberOfPoints())
{
vtkErrorMacro(<< "point id " << id << " is larger than the number of points");
return;
}
this->Points->SetPoint(id, x, y, z);
this->Modified();
}
//----------------------------------------------------------------------------
void vtkPolyLineSource::SetPoints(vtkPoints* points)
{
if ( points != this->Points )
{
if ( this->Points != nullptr )
{
this->Points->Delete();
}
this->Points = points;
if ( this->Points != nullptr )
{
this->Points->Register(this);
}
this->Modified();
}
}
//----------------------------------------------------------------------------
......@@ -161,6 +74,5 @@ void vtkPolyLineSource::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
os << indent << "Points: " << this->Points << "\n";
os << indent << "Closed: " << this->Closed << "\n";
}
......@@ -24,43 +24,17 @@
#define vtkPolyLineSource_h
#include "vtkFiltersSourcesModule.h" // For export macro
#include "vtkPolyDataAlgorithm.h"
#include "vtkPolyPointSource.h"
class vtkPoints;
class VTKFILTERSSOURCES_EXPORT vtkPolyLineSource : public vtkPolyDataAlgorithm
class VTKFILTERSSOURCES_EXPORT vtkPolyLineSource : public vtkPolyPointSource
{
public:
static vtkPolyLineSource* New();
vtkTypeMacro(vtkPolyLineSource, vtkPolyDataAlgorithm);
vtkTypeMacro(vtkPolyLineSource, vtkPolyPointSource);
void PrintSelf(ostream& os, vtkIndent indent) override;
//@{
/**
* Set the number of points in the poly line.
*/
void SetNumberOfPoints(vtkIdType numPoints);
vtkIdType GetNumberOfPoints();
//@}
/**
* Resize while preserving data.
*/
void Resize(vtkIdType numPoints);
/**
* Set a point location.
*/
void SetPoint(vtkIdType id, double x, double y, double z);
//@{
/**
* Get the points.
*/
void SetPoints(vtkPoints* points);
vtkGetObjectMacro(Points, vtkPoints);
//@}
//@{
/**
* Set whether to close the poly line by connecting the last and first points.
......@@ -76,8 +50,6 @@ protected:
int RequestData(vtkInformation*, vtkInformationVector**, vtkInformationVector *) override;
vtkPoints* Points;
vtkTypeBool Closed;
private:
......
/*=========================================================================
Program: Visualization Toolkit
Module: vtkPolyPointSource.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 "vtkPolyPointSource.h"
#include "vtkCellArray.h"
#include "vtkIdList.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkObjectFactory.h"
#include "vtkSmartPointer.h"
//----------------------------------------------------------------------------
vtkStandardNewMacro(vtkPolyPointSource);
vtkCxxSetObjectMacro(vtkPolyPointSource, Points, vtkPoints)
//----------------------------------------------------------------------------
vtkPolyPointSource::vtkPolyPointSource()
{
this->Points = nullptr;
this->SetNumberOfInputPorts(0);
}
//----------------------------------------------------------------------------
vtkPolyPointSource::~vtkPolyPointSource()
{
if (this->Points)
{
this->Points->Delete();
}
}
//----------------------------------------------------------------------------
vtkMTimeType vtkPolyPointSource::GetMTime()
{
vtkMTimeType mTime = this->Superclass::GetMTime();
vtkMTimeType time;
if ( this->Points != nullptr )
{
time = this->Points->GetMTime();
mTime = ( time > mTime ? time : mTime );
}
return mTime;
}
//----------------------------------------------------------------------------
void vtkPolyPointSource::SetNumberOfPoints(vtkIdType numPoints)
{
if (!this->Points)
{
vtkPoints* pts = vtkPoints::New(VTK_DOUBLE);
this->SetPoints(pts);
this->Points = pts;
pts->Delete();
}
if (numPoints != this->GetNumberOfPoints())
{
this->Points->SetNumberOfPoints(numPoints);
this->Modified();
}
}
//----------------------------------------------------------------------------
vtkIdType vtkPolyPointSource::GetNumberOfPoints()
{
if (this->Points)
{
return this->Points->GetNumberOfPoints();
}
return 0;
}
//----------------------------------------------------------------------------
void vtkPolyPointSource::Resize(vtkIdType numPoints)
{
if (!this->Points)
{
this->SetNumberOfPoints(numPoints);
}
if (numPoints != this->GetNumberOfPoints())
{
this->Points->Resize(numPoints);
this->Modified();
}
}
//----------------------------------------------------------------------------
void vtkPolyPointSource::SetPoint(vtkIdType id, double x, double y, double z)
{
if (!this->Points)
{
return;
}
if (id >= this->Points->GetNumberOfPoints())
{
vtkErrorMacro(<< "point id " << id << " is larger than the number of points");
return;
}
this->Points->SetPoint(id, x, y, z);
this->Modified();
}
//----------------------------------------------------------------------------
int vtkPolyPointSource::RequestData(
vtkInformation *vtkNotUsed(request),
vtkInformationVector **vtkNotUsed(inputVector),
vtkInformationVector *outputVector)
{
// get the info object
vtkInformation *outInfo = outputVector->GetInformationObject(0);
// get the output
vtkPolyData *output = vtkPolyData::SafeDownCast(
outInfo->Get(vtkDataObject::DATA_OBJECT()));
vtkIdType numPoints = this->GetNumberOfPoints();
vtkSmartPointer<vtkIdList> pointIds = vtkSmartPointer<vtkIdList>::New();
pointIds->SetNumberOfIds(numPoints);
for (vtkIdType i = 0; i < numPoints; ++i)
{
pointIds->SetId(i, i);
}
vtkSmartPointer<vtkCellArray> PolyPoint = vtkSmartPointer<vtkCellArray>::New();
PolyPoint->InsertNextCell(pointIds);
output->SetPoints(this->Points);
output->SetLines(PolyPoint);
return 1;
}
//----------------------------------------------------------------------------
void vtkPolyPointSource::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
os << indent << "Points: " << this->Points << "\n";
}
/*=========================================================================
Program: Visualization Toolkit
Module: vtkPolyPointSource.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.
=========================================================================*/
/**
* @class vtkPolyPointSource
* @brief create points from a list of input points
*
* vtkPolyPointSource is a source object that creates a vert from
* user-specified points. The output is a vtkPolyData.
*/
#ifndef vtkPolyPointSource_h
#define vtkPolyPointSource_h
#include "vtkFiltersSourcesModule.h" // For export macro
#include "vtkPolyDataAlgorithm.h"
class vtkPoints;
class VTKFILTERSSOURCES_EXPORT vtkPolyPointSource : public vtkPolyDataAlgorithm
{
public:
static vtkPolyPointSource* New();
vtkTypeMacro(vtkPolyPointSource, vtkPolyDataAlgorithm);
void PrintSelf(ostream& os, vtkIndent indent) override;
//@{
/**
* Set the number of points in the poly line.
*/
void SetNumberOfPoints(vtkIdType numPoints);
vtkIdType GetNumberOfPoints();
//@}
/**
* Resize while preserving data.
*/
void Resize(vtkIdType numPoints);
/**
* Set a point location.
*/
void SetPoint(vtkIdType id, double x, double y, double z);
//@{
/**
* Get the points.
*/
void SetPoints(vtkPoints* points);
vtkGetObjectMacro(Points, vtkPoints);
//@}
/**
* Get the mtime plus consider its Points
*/
vtkMTimeType GetMTime() override;
protected:
vtkPolyPointSource();
~vtkPolyPointSource() override;
int RequestData(vtkInformation*, vtkInformationVector**, vtkInformationVector *) override;
vtkPoints* Points;
vtkTypeBool Closed;
private:
vtkPolyPointSource(const vtkPolyPointSource&) = delete;
void operator=(const vtkPolyPointSource&) = delete;
};
#endif
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