OGS
VtkPolylinesSource.cpp
Go to the documentation of this file.
1
15// ** INCLUDES **
16#include "VtkPolylinesSource.h"
17
18#include <vtkCellArray.h>
19#include <vtkCellData.h>
20#include <vtkInformation.h>
21#include <vtkInformationVector.h>
22#include <vtkObjectFactory.h>
23#include <vtkPointData.h>
24#include <vtkPoints.h>
25#include <vtkPolyData.h>
26#include <vtkProperty.h>
27#include <vtkSmartPointer.h>
28#include <vtkStreamingDemandDrivenPipeline.h>
29
31#include "BaseLib/Logging.h"
32#include "GeoLib/Polyline.h"
33
35
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}
44
46
47void VtkPolylinesSource::PrintSelf(ostream& os, vtkIndent indent)
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}
69
70int VtkPolylinesSource::RequestData(vtkInformation* request,
71 vtkInformationVector** inputVector,
72 vtkInformationVector* outputVector)
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}
150
152 vtkInformation* /*request*/,
153 vtkInformationVector** /*inputVector*/,
154 vtkInformationVector* /*outputVector*/)
155{
156 return 1;
157}
158
159void VtkPolylinesSource::SetUserProperty(QString name, QVariant value)
160{
161 Q_UNUSED(name);
162 Q_UNUSED(value);
163}
Definition of the Color class.
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
Definition of the PolyLine class.
vtkStandardNewMacro(VtkPolylinesSource)
Definition of the VtkPolylinesSource class.
const double * data() const
Definition Point3d.h:60
vtkProperty * GetProperties() const
Returns the vtk properties.
VtkPolylinesSource is a VTK source object for the visualisation of polyline data. As a vtkPolyDataAlg...
int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
~VtkPolylinesSource() override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Computes the polygonal data object.
const std::vector< GeoLib::Polyline * > * _polylines
The polylines to visualize.
void SetUserProperty(QString name, QVariant value) override
Sets a user property. This should be implemented by subclasses.
void PrintSelf(ostream &os, vtkIndent indent) override
Prints its data on a stream.
Color getRandomColor()
Returns a random RGB colour.
Definition Color.cpp:29
std::array< unsigned char, 4 > Color
Definition Color.h:24