OGS
VtkPolylinesSource Class Reference

Detailed Description

VtkPolylinesSource is a VTK source object for the visualisation of polyline data. As a vtkPolyDataAlgorithm it outputs polygonal data.

Definition at line 20 of file VtkPolylinesSource.h.

#include <VtkPolylinesSource.h>

Inheritance diagram for VtkPolylinesSource:
[legend]
Collaboration diagram for VtkPolylinesSource:
[legend]

Public Member Functions

 vtkTypeMacro (VtkPolylinesSource, vtkPolyDataAlgorithm)
void setPolylines (const std::vector< GeoLib::Polyline * > *polylines)
 Sets the polyline vector.
void PrintSelf (ostream &os, vtkIndent indent) override
 Prints its data on a stream.
void SetUserProperty (QString name, QVariant value) override
 Sets a user property. This should be implemented by subclasses.
Public Member Functions inherited from VtkAlgorithmProperties
 VtkAlgorithmProperties (QObject *parent=nullptr)
 Constructor (sets default values)
 ~VtkAlgorithmProperties () override
vtkProperty * GetProperties () const
 Returns the vtk properties.
vtkTexture * GetTexture ()
 Returns a texture (if one has been assigned).
void SetTexture (vtkTexture *t)
 Sets a texture for the VtkVisPipelineItem.
vtkLookupTable * GetLookupTable (const QString &array_name)
 Returns the colour lookup table (if one has been assigned).
void RemoveLookupTable (const QString &array_name)
 Removes the lookup table for the given scalar.
void SetLookUpTable (const QString &array_name, vtkLookupTable *lut)
 Sets a colour lookup table for the given scalar array of the VtkVisPipelineItem.
void SetLookUpTable (const QString &array_name, const QString &filename)
 Loads a predefined color lookup table from a file for the specified scalar array.
bool GetScalarVisibility () const
 Returns the scalar visibility.
void SetScalarVisibility (bool on)
 Sets the scalar visibility.
QString GetName () const
 Returns the name. This is set to the file path if it is a source algorithm.
void SetName (QString name)
 Sets the name.
bool IsRemovable () const
 Is this algorithm removable from the pipeline (view).
QMap< QString, QVariant > * GetAlgorithmUserProperties () const
 Returns a map of user properties.
QMap< QString, QList< QVariant > > * GetAlgorithmUserVectorProperties () const
 Returns a map of vector user properties.
QVariant GetUserProperty (QString name) const
 Returns the value of a user property.
virtual void SetUserVectorProperty (QString name, QList< QVariant > values)
 Sets a vector user property. This should be implemented by subclasses.
QList< QVariant > GetUserVectorProperty (QString name) const
 Returns a list of values of a vector user property.
void SetActiveAttribute (QString name)
 Set the active attribute.
QString GetActiveAttribute () const
 Returns the desired active attribute.

Static Public Member Functions

static VtkPolylinesSourceNew ()
 Create new objects with New() because of VTKs object reference counting.

Protected Member Functions

 VtkPolylinesSource ()
 ~VtkPolylinesSource () override
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 Computes the polygonal data object.
int RequestInformation (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override

Protected Attributes

const std::vector< GeoLib::Polyline * > * _polylines {nullptr}
 The polylines to visualize.
Protected Attributes inherited from VtkAlgorithmProperties
vtkProperty * _property
vtkTexture * _texture
bool _scalarVisibility
std::map< QString, vtkLookupTable * > _lut
QString _name
QString _activeAttributeName
bool _removable
QMap< QString, QVariant > * _algorithmUserProperties
QMap< QString, QList< QVariant > > * _algorithmUserVectorProperties

Additional Inherited Members

Signals inherited from VtkAlgorithmProperties
void ScalarVisibilityChanged (bool on)

Constructor & Destructor Documentation

◆ VtkPolylinesSource()

VtkPolylinesSource::VtkPolylinesSource ( )
protected

Definition at line 25 of file VtkPolylinesSource.cpp.

26{
27 _removable = false; // From VtkAlgorithmProperties
28 this->SetNumberOfInputPorts(0);
29
31 GetProperties()->SetColor(c[0] / 255.0, c[1] / 255.0, c[2] / 255.0);
32}
vtkProperty * GetProperties() const
Returns the vtk properties.
Color getRandomColor()
Returns a random RGB colour.
Definition Color.cpp:18
std::array< unsigned char, 4 > Color
Definition Color.h:12

References VtkAlgorithmProperties::_removable, VtkAlgorithmProperties::GetProperties(), and DataHolderLib::getRandomColor().

Referenced by New(), and vtkTypeMacro().

◆ ~VtkPolylinesSource()

VtkPolylinesSource::~VtkPolylinesSource ( )
overrideprotecteddefault

Member Function Documentation

◆ New()

VtkPolylinesSource * VtkPolylinesSource::New ( )
static

Create new objects with New() because of VTKs object reference counting.

References VtkPolylinesSource().

◆ PrintSelf()

void VtkPolylinesSource::PrintSelf ( ostream & os,
vtkIndent indent )
override

Prints its data on a stream.

Definition at line 36 of file VtkPolylinesSource.cpp.

37{
38 this->Superclass::PrintSelf(os, indent);
39
40 if (_polylines->empty())
41 {
42 return;
43 }
44
45 for (auto polyline : *_polylines)
46 {
47 os << indent << "== Polyline =="
48 << "\n";
49 int numPoints = polyline->getNumberOfPoints();
50 for (int i = 0; i < numPoints; i++)
51 {
52 const GeoLib::Point& point = *polyline->getPoint(i);
53 os << indent << "Point " << i << " (" << point[0] << ", "
54 << point[1] << ", " << point[2] << ")\n";
55 }
56 }
57}
const std::vector< GeoLib::Polyline * > * _polylines
The polylines to visualize.

References _polylines.

◆ RequestData()

int VtkPolylinesSource::RequestData ( vtkInformation * request,
vtkInformationVector ** inputVector,
vtkInformationVector * outputVector )
overrideprotected

Computes the polygonal data object.

Definition at line 59 of file VtkPolylinesSource.cpp.

62{
63 (void)request;
64 (void)inputVector;
65
66 if (!_polylines)
67 {
68 return 0;
69 }
70 if (_polylines->empty())
71 {
72 ERR("VtkPolylineSource::RequestData(): Size of polyline vector is 0");
73 return 0;
74 }
75
76 vtkSmartPointer<vtkInformation> outInfo =
77 outputVector->GetInformationObject(0);
78 vtkSmartPointer<vtkPolyData> output =
79 vtkPolyData::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
80
81 vtkSmartPointer<vtkPoints> newPoints = vtkSmartPointer<vtkPoints>::New();
82 vtkSmartPointer<vtkCellArray> newLines =
83 vtkSmartPointer<vtkCellArray>::New();
84
85 // newPoints->Allocate(numPoints);
86 // newLines->Allocate(newLines->EstimateSize(numLines, 2));
87
88 if (outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER()) >
89 0)
90 {
91 return 1;
92 }
93
94 vtkSmartPointer<vtkIntArray> plyIDs = vtkSmartPointer<vtkIntArray>::New();
95 plyIDs->SetNumberOfComponents(1);
96 plyIDs->SetName("PolylineIDs");
97
98 unsigned lastMaxIndex(0);
99 const std::size_t nPolylines(_polylines->size());
100 for (std::size_t j = 0; j < nPolylines; j++)
101 {
102 const int numPoints = (*_polylines)[j]->getNumberOfPoints();
103 const bool isClosed = (*_polylines)[j]->isClosed();
104
105 // Generate points
106 const int numVerts = (isClosed) ? numPoints - 1 : numPoints;
107 for (int i = 0; i < numVerts; i++)
108 {
109 const GeoLib::Point* point = (*_polylines)[j]->getPoint(i);
110 const double* coords = point->data();
111 newPoints->InsertNextPoint(coords);
112 }
113
114 // Generate lines
115 newLines->InsertNextCell(numPoints);
116 plyIDs->InsertNextValue(j);
117 for (int i = 0; i < numVerts; i++)
118 {
119 newLines->InsertCellPoint(i + lastMaxIndex);
120 }
121
122 if (isClosed)
123 {
124 newLines->InsertCellPoint(lastMaxIndex);
125 }
126
127 lastMaxIndex += numVerts;
128 }
129
130 output->SetPoints(newPoints);
131 output->SetLines(newLines);
132 output->GetCellData()->AddArray(plyIDs);
133 output->GetCellData()->SetActiveAttribute("PolylineIDs",
134 vtkDataSetAttributes::SCALARS);
135 output->Squeeze();
136
137 return 1;
138}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
const double * data() const
Definition Point3d.h:51
constexpr ranges::views::view_closure coords
Definition Mesh.h:223

References _polylines, MathLib::Point3d::data(), and ERR().

◆ RequestInformation()

int VtkPolylinesSource::RequestInformation ( vtkInformation * request,
vtkInformationVector ** inputVector,
vtkInformationVector * outputVector )
overrideprotected

Definition at line 140 of file VtkPolylinesSource.cpp.

144{
145 return 1;
146}

◆ setPolylines()

void VtkPolylinesSource::setPolylines ( const std::vector< GeoLib::Polyline * > * polylines)
inline

Sets the polyline vector.

Definition at line 29 of file VtkPolylinesSource.h.

29{ _polylines = polylines; }

References _polylines.

Referenced by GeoObjectListItem::GeoObjectListItem().

◆ SetUserProperty()

void VtkPolylinesSource::SetUserProperty ( QString name,
QVariant value )
overridevirtual

Sets a user property. This should be implemented by subclasses.

Reimplemented from VtkAlgorithmProperties.

Definition at line 148 of file VtkPolylinesSource.cpp.

149{
150 Q_UNUSED(name);
151 Q_UNUSED(value);
152}

◆ vtkTypeMacro()

VtkPolylinesSource::vtkTypeMacro ( VtkPolylinesSource ,
vtkPolyDataAlgorithm  )

References VtkPolylinesSource().

Member Data Documentation

◆ _polylines

const std::vector<GeoLib::Polyline*>* VtkPolylinesSource::_polylines {nullptr}
protected

The polylines to visualize.

Definition at line 50 of file VtkPolylinesSource.h.

50{nullptr};

Referenced by PrintSelf(), RequestData(), and setPolylines().


The documentation for this class was generated from the following files: