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 31 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 36 of file VtkPolylinesSource.cpp.

37{
38 _removable = false; // From VtkAlgorithmProperties
39 this->SetNumberOfInputPorts(0);
40
42 GetProperties()->SetColor(c[0] / 255.0, c[1] / 255.0, c[2] / 255.0);
43}
vtkProperty * GetProperties() const
Returns the vtk properties.
Color getRandomColor()
Returns a random RGB colour.
Definition Color.cpp:29
std::array< unsigned char, 4 > Color
Definition Color.h:24

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

◆ ~VtkPolylinesSource()

VtkPolylinesSource::~VtkPolylinesSource ( )
overrideprotecteddefault

Member Function Documentation

◆ New()

static VtkPolylinesSource * VtkPolylinesSource::New ( )
static

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

◆ PrintSelf()

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

Prints its data on a stream.

Definition at line 47 of file VtkPolylinesSource.cpp.

48{
49 this->Superclass::PrintSelf(os, indent);
50
51 if (_polylines->empty())
52 {
53 return;
54 }
55
56 for (auto polyline : *_polylines)
57 {
58 os << indent << "== Polyline =="
59 << "\n";
60 int numPoints = polyline->getNumberOfPoints();
61 for (int i = 0; i < numPoints; i++)
62 {
63 const GeoLib::Point& point = *polyline->getPoint(i);
64 os << indent << "Point " << i << " (" << point[0] << ", "
65 << point[1] << ", " << point[2] << ")\n";
66 }
67 }
68}
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 70 of file VtkPolylinesSource.cpp.

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

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

◆ RequestInformation()

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

Definition at line 151 of file VtkPolylinesSource.cpp.

155{
156 return 1;
157}

◆ setPolylines()

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

Sets the polyline vector.

Definition at line 40 of file VtkPolylinesSource.h.

40{ _polylines = polylines; }

References _polylines.

◆ 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 159 of file VtkPolylinesSource.cpp.

160{
161 Q_UNUSED(name);
162 Q_UNUSED(value);
163}

◆ vtkTypeMacro()

VtkPolylinesSource::vtkTypeMacro ( VtkPolylinesSource ,
vtkPolyDataAlgorithm  )

Member Data Documentation

◆ _polylines

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

The polylines to visualize.

Definition at line 61 of file VtkPolylinesSource.h.

61{nullptr};

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


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