OGS
VtkSurfacesSource.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 **
5#include "VtkSurfacesSource.h"
6
7#include <vtkCellArray.h>
8#include <vtkCellData.h>
9#include <vtkInformation.h>
10#include <vtkInformationVector.h>
11#include <vtkObjectFactory.h>
12#include <vtkPolyData.h>
13#include <vtkProperty.h>
14#include <vtkSmartPointer.h>
15#include <vtkStreamingDemandDrivenPipeline.h>
16#include <vtkTriangle.h>
17
18#include <limits>
19
21#include "GeoLib/Triangle.h"
22
24
26{
27 _removable = false; // From VtkAlgorithmProperties
28 this->SetNumberOfInputPorts(0);
29 // this->SetColorBySurface(true);
30
32 vtkProperty* vtkProps = GetProperties();
33 vtkProps->SetColor(c[0] / 255.0, c[1] / 255.0, c[2] / 255.0);
34 vtkProps->SetEdgeVisibility(0);
35}
36
37void VtkSurfacesSource::PrintSelf(ostream& os, vtkIndent indent)
38{
39 this->Superclass::PrintSelf(os, indent);
40
41 if (_surfaces->empty())
42 {
43 return;
44 }
45
46 os << indent << "== VtkSurfacesSource =="
47 << "\n";
48}
49
50int VtkSurfacesSource::RequestData(vtkInformation* request,
51 vtkInformationVector** inputVector,
52 vtkInformationVector* outputVector)
53{
54 (void)request;
55 (void)inputVector;
56
57 if (_surfaces->empty())
58 {
59 return 0;
60 }
61
62 const std::vector<GeoLib::Point*>* surfacePoints =
63 (*_surfaces)[0]->getPointVec();
64 std::size_t const nPoints = surfacePoints->size();
65
66 vtkSmartPointer<vtkInformation> outInfo =
67 outputVector->GetInformationObject(0);
68 vtkSmartPointer<vtkPolyData> output =
69 vtkPolyData::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
70 if (outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER()) >
71 0)
72 {
73 return 1;
74 }
75
76 vtkSmartPointer<vtkPoints> newPoints = vtkSmartPointer<vtkPoints>::New();
77 newPoints->SetNumberOfPoints(nPoints);
78
79 vtkSmartPointer<vtkCellArray> newPolygons =
80 vtkSmartPointer<vtkCellArray>::New();
81 // newPolygons->Allocate(nSurfaces);
82
83 vtkSmartPointer<vtkIntArray> sfcIDs = vtkSmartPointer<vtkIntArray>::New();
84 sfcIDs->SetNumberOfComponents(1);
85 sfcIDs->SetName("SurfaceIDs");
86
87 for (std::size_t i = 0; i < nPoints; ++i)
88 {
89 newPoints->SetPoint(i, (*surfacePoints)[i]->data());
90 }
91
92 int count(0);
93 for (auto surface : *_surfaces)
94 {
95 const std::size_t nTriangles = surface->getNumberOfTriangles();
96
97 for (std::size_t i = 0; i < nTriangles; ++i)
98 {
99 vtkTriangle* new_tri = vtkTriangle::New();
100 new_tri->GetPointIds()->SetNumberOfIds(3);
101
102 const GeoLib::Triangle* triangle = (*surface)[i];
103 for (std::size_t j = 0; j < 3; ++j)
104 {
105 new_tri->GetPointIds()->SetId(j, ((*triangle)[j]));
106 }
107 newPolygons->InsertNextCell(new_tri);
108 sfcIDs->InsertNextValue(count);
109 new_tri->Delete();
110 }
111 count++;
112 }
113
114 output->SetPoints(newPoints);
115 output->SetPolys(newPolygons);
116 output->GetCellData()->AddArray(sfcIDs);
117 output->GetCellData()->SetActiveAttribute("SurfaceIDs",
118 vtkDataSetAttributes::SCALARS);
119 output->Squeeze();
120
121 return 1;
122}
123
125 vtkInformation* /*request*/,
126 vtkInformationVector** /*inputVector*/,
127 vtkInformationVector* /*outputVector*/)
128{
129 return 1;
130}
131
132void VtkSurfacesSource::SetUserProperty(QString name, QVariant value)
133{
135 (*_algorithmUserProperties)[name] = value;
136}
vtkStandardNewMacro(VtkSurfacesSource)
Class Triangle consists of a reference to a point vector and a vector that stores the indices in the ...
Definition Triangle.h:21
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:18
std::array< unsigned char, 4 > Color
Definition Color.h:12