OGS
VtkPolylinesSource.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4// ** INCLUDES **
6
7#include <vtkCellArray.h>
8#include <vtkCellData.h>
9#include <vtkInformation.h>
10#include <vtkInformationVector.h>
11#include <vtkObjectFactory.h>
12#include <vtkPointData.h>
13#include <vtkPoints.h>
14#include <vtkPolyData.h>
15#include <vtkProperty.h>
16#include <vtkSmartPointer.h>
17#include <vtkStreamingDemandDrivenPipeline.h>
18
20#include "BaseLib/Logging.h"
21#include "GeoLib/Polyline.h"
22
24
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}
33
35
36void VtkPolylinesSource::PrintSelf(ostream& os, vtkIndent indent)
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}
58
59int VtkPolylinesSource::RequestData(vtkInformation* request,
60 vtkInformationVector** inputVector,
61 vtkInformationVector* outputVector)
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}
139
141 vtkInformation* /*request*/,
142 vtkInformationVector** /*inputVector*/,
143 vtkInformationVector* /*outputVector*/)
144{
145 return 1;
146}
147
148void VtkPolylinesSource::SetUserProperty(QString name, QVariant value)
149{
150 Q_UNUSED(name);
151 Q_UNUSED(value);
152}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
vtkStandardNewMacro(VtkPolylinesSource)
const double * data() const
Definition Point3d.h:51
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:18
std::array< unsigned char, 4 > Color
Definition Color.h:12