OGS
VtkSurfacesSource.cpp
Go to the documentation of this file.
1
15// ** INCLUDES **
16#include "VtkSurfacesSource.h"
17
18#include <vtkCellArray.h>
19#include <vtkCellData.h>
20#include <vtkInformation.h>
21#include <vtkInformationVector.h>
22#include <vtkObjectFactory.h>
23#include <vtkPolyData.h>
24#include <vtkProperty.h>
25#include <vtkSmartPointer.h>
26#include <vtkStreamingDemandDrivenPipeline.h>
27#include <vtkTriangle.h>
28
29#include <limits>
30
32#include "GeoLib/Triangle.h"
33
35
37{
38 _removable = false; // From VtkAlgorithmProperties
39 this->SetNumberOfInputPorts(0);
40 // this->SetColorBySurface(true);
41
43 vtkProperty* vtkProps = GetProperties();
44 vtkProps->SetColor(c[0] / 255.0, c[1] / 255.0, c[2] / 255.0);
45 vtkProps->SetEdgeVisibility(0);
46}
47
48void VtkSurfacesSource::PrintSelf(ostream& os, vtkIndent indent)
49{
50 this->Superclass::PrintSelf(os, indent);
51
52 if (_surfaces->empty())
53 {
54 return;
55 }
56
57 os << indent << "== VtkSurfacesSource =="
58 << "\n";
59}
60
61int VtkSurfacesSource::RequestData(vtkInformation* request,
62 vtkInformationVector** inputVector,
63 vtkInformationVector* outputVector)
64{
65 (void)request;
66 (void)inputVector;
67
68 if (_surfaces->empty())
69 {
70 return 0;
71 }
72
73 const std::vector<GeoLib::Point*>* surfacePoints =
74 (*_surfaces)[0]->getPointVec();
75 std::size_t const nPoints = surfacePoints->size();
76
77 vtkSmartPointer<vtkInformation> outInfo =
78 outputVector->GetInformationObject(0);
79 vtkSmartPointer<vtkPolyData> output =
80 vtkPolyData::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
81 if (outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER()) >
82 0)
83 {
84 return 1;
85 }
86
87 vtkSmartPointer<vtkPoints> newPoints = vtkSmartPointer<vtkPoints>::New();
88 newPoints->SetNumberOfPoints(nPoints);
89
90 vtkSmartPointer<vtkCellArray> newPolygons =
91 vtkSmartPointer<vtkCellArray>::New();
92 // newPolygons->Allocate(nSurfaces);
93
94 vtkSmartPointer<vtkIntArray> sfcIDs = vtkSmartPointer<vtkIntArray>::New();
95 sfcIDs->SetNumberOfComponents(1);
96 sfcIDs->SetName("SurfaceIDs");
97
98 for (std::size_t i = 0; i < nPoints; ++i)
99 {
100 newPoints->SetPoint(i, (*surfacePoints)[i]->data());
101 }
102
103 int count(0);
104 for (auto surface : *_surfaces)
105 {
106 const std::size_t nTriangles = surface->getNumberOfTriangles();
107
108 for (std::size_t i = 0; i < nTriangles; ++i)
109 {
110 vtkTriangle* new_tri = vtkTriangle::New();
111 new_tri->GetPointIds()->SetNumberOfIds(3);
112
113 const GeoLib::Triangle* triangle = (*surface)[i];
114 for (std::size_t j = 0; j < 3; ++j)
115 {
116 new_tri->GetPointIds()->SetId(j, ((*triangle)[j]));
117 }
118 newPolygons->InsertNextCell(new_tri);
119 sfcIDs->InsertNextValue(count);
120 new_tri->Delete();
121 }
122 count++;
123 }
124
125 output->SetPoints(newPoints);
126 output->SetPolys(newPolygons);
127 output->GetCellData()->AddArray(sfcIDs);
128 output->GetCellData()->SetActiveAttribute("SurfaceIDs",
129 vtkDataSetAttributes::SCALARS);
130 output->Squeeze();
131
132 return 1;
133}
134
136 vtkInformation* /*request*/,
137 vtkInformationVector** /*inputVector*/,
138 vtkInformationVector* /*outputVector*/)
139{
140 return 1;
141}
142
143void VtkSurfacesSource::SetUserProperty(QString name, QVariant value)
144{
146 (*_algorithmUserProperties)[name] = value;
147}
Definition of the Color class.
vtkStandardNewMacro(VtkSurfacesSource)
Definition of the VtkSurfacesSource class.
Class Triangle consists of a reference to a point vector and a vector that stores the indices in the ...
Definition Triangle.h:27
vtkProperty * GetProperties() const
Returns the vtk properties.
virtual void SetUserProperty(QString name, QVariant value)
Sets a user property. This should be implemented by subclasses.
VTK source object for the visualisation of surfaces. Technically, surfaces are displayed as triangula...
int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
void PrintSelf(ostream &os, vtkIndent indent) override
Prints its data on a stream.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Computes the polygonal data object.
const std::vector< GeoLib::Surface * > * _surfaces
The surfaces to visualize.
void SetUserProperty(QString name, QVariant value) override
Generates random colors for each surface.
Color getRandomColor()
Returns a random RGB colour.
Definition Color.cpp:29
std::array< unsigned char, 4 > Color
Definition Color.h:24