Loading [MathJax]/extensions/tex2jax.js
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 
48 void 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 
61 int VtkSurfacesSource::RequestData(vtkInformation* request,
62  vtkInformationVector** inputVector,
63  vtkInformationVector* outputVector)
64 {
65  (void)request;
66  (void)inputVector;
67 
68  const int nSurfaces = _surfaces->size();
69  if (nSurfaces == 0)
70  {
71  return 0;
72  }
73 
74  const std::vector<GeoLib::Point*>* surfacePoints =
75  (*_surfaces)[0]->getPointVec();
76  std::size_t nPoints = surfacePoints->size();
77 
78  vtkSmartPointer<vtkInformation> outInfo =
79  outputVector->GetInformationObject(0);
80  vtkSmartPointer<vtkPolyData> output =
81  vtkPolyData::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
82  if (outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER()) >
83  0)
84  {
85  return 1;
86  }
87 
88  vtkSmartPointer<vtkPoints> newPoints = vtkSmartPointer<vtkPoints>::New();
89  newPoints->SetNumberOfPoints(nPoints);
90 
91  vtkSmartPointer<vtkCellArray> newPolygons =
92  vtkSmartPointer<vtkCellArray>::New();
93  // newPolygons->Allocate(nSurfaces);
94 
95  vtkSmartPointer<vtkIntArray> sfcIDs = vtkSmartPointer<vtkIntArray>::New();
96  sfcIDs->SetNumberOfComponents(1);
97  sfcIDs->SetName("SurfaceIDs");
98 
99  for (std::size_t i = 0; i < nPoints; ++i)
100  {
101  const double* coords =
102  const_cast<double*>((*surfacePoints)[i]->getCoords());
103  newPoints->SetPoint(i, coords);
104  }
105 
106  int count(0);
107  for (auto surface : *_surfaces)
108  {
109  const std::size_t nTriangles = surface->getNumberOfTriangles();
110 
111  for (std::size_t i = 0; i < nTriangles; ++i)
112  {
113  vtkTriangle* new_tri = vtkTriangle::New();
114  new_tri->GetPointIds()->SetNumberOfIds(3);
115 
116  const GeoLib::Triangle* triangle = (*surface)[i];
117  for (std::size_t j = 0; j < 3; ++j)
118  {
119  new_tri->GetPointIds()->SetId(j, ((*triangle)[j]));
120  }
121  newPolygons->InsertNextCell(new_tri);
122  sfcIDs->InsertNextValue(count);
123  new_tri->Delete();
124  }
125  count++;
126  }
127 
128  output->SetPoints(newPoints);
129  output->SetPolys(newPolygons);
130  output->GetCellData()->AddArray(sfcIDs);
131  output->GetCellData()->SetActiveAttribute("SurfaceIDs",
132  vtkDataSetAttributes::SCALARS);
133  output->Squeeze();
134 
135  return 1;
136 }
137 
139  vtkInformation* /*request*/,
140  vtkInformationVector** /*inputVector*/,
141  vtkInformationVector* /*outputVector*/)
142 {
143  return 1;
144 }
145 
146 void VtkSurfacesSource::SetUserProperty(QString name, QVariant value)
147 {
149  (*_algorithmUserProperties)[name] = value;
150 }
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:26
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.
std::array< unsigned char, 4 > Color
Definition: Color.h:24
Color getRandomColor()
Returns a random RGB colour.
Definition: Color.cpp:29