OGS
VtkSurfacesSource Class Reference

Detailed Description

VTK source object for the visualisation of surfaces. Technically, surfaces are displayed as triangulated polydata objects.

Definition at line 16 of file VtkSurfacesSource.h.

#include <VtkSurfacesSource.h>

Inheritance diagram for VtkSurfacesSource:
[legend]
Collaboration diagram for VtkSurfacesSource:
[legend]

Public Member Functions

 vtkTypeMacro (VtkSurfacesSource, vtkPolyDataAlgorithm)
void setSurfaces (const std::vector< GeoLib::Surface * > *surfaces)
 Sets the surfaces vector.
void PrintSelf (ostream &os, vtkIndent indent) override
 Prints its data on a stream.
void SetUserProperty (QString name, QVariant value) override
 Generates random colors for each surface.
Public Member Functions inherited from VtkAlgorithmProperties
 VtkAlgorithmProperties (QObject *parent=nullptr)
 Constructor (sets default values)
 ~VtkAlgorithmProperties () override
vtkProperty * GetProperties () const
 Returns the vtk properties.
vtkTexture * GetTexture ()
 Returns a texture (if one has been assigned).
void SetTexture (vtkTexture *t)
 Sets a texture for the VtkVisPipelineItem.
vtkLookupTable * GetLookupTable (const QString &array_name)
 Returns the colour lookup table (if one has been assigned).
void RemoveLookupTable (const QString &array_name)
 Removes the lookup table for the given scalar.
void SetLookUpTable (const QString &array_name, vtkLookupTable *lut)
 Sets a colour lookup table for the given scalar array of the VtkVisPipelineItem.
void SetLookUpTable (const QString &array_name, const QString &filename)
 Loads a predefined color lookup table from a file for the specified scalar array.
bool GetScalarVisibility () const
 Returns the scalar visibility.
void SetScalarVisibility (bool on)
 Sets the scalar visibility.
QString GetName () const
 Returns the name. This is set to the file path if it is a source algorithm.
void SetName (QString name)
 Sets the name.
bool IsRemovable () const
 Is this algorithm removable from the pipeline (view).
QMap< QString, QVariant > * GetAlgorithmUserProperties () const
 Returns a map of user properties.
QMap< QString, QList< QVariant > > * GetAlgorithmUserVectorProperties () const
 Returns a map of vector user properties.
QVariant GetUserProperty (QString name) const
 Returns the value of a user property.
virtual void SetUserVectorProperty (QString name, QList< QVariant > values)
 Sets a vector user property. This should be implemented by subclasses.
QList< QVariant > GetUserVectorProperty (QString name) const
 Returns a list of values of a vector user property.
void SetActiveAttribute (QString name)
 Set the active attribute.
QString GetActiveAttribute () const
 Returns the desired active attribute.

Static Public Member Functions

static VtkSurfacesSourceNew ()
 Create new objects with New() because of VTKs object reference counting.

Protected Member Functions

 VtkSurfacesSource ()
 ~VtkSurfacesSource () override=default
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 Computes the polygonal data object.
int RequestInformation (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override

Protected Attributes

const std::vector< GeoLib::Surface * > * _surfaces {nullptr}
 The surfaces to visualize.
Protected Attributes inherited from VtkAlgorithmProperties
vtkProperty * _property
vtkTexture * _texture
bool _scalarVisibility
std::map< QString, vtkLookupTable * > _lut
QString _name
QString _activeAttributeName
bool _removable
QMap< QString, QVariant > * _algorithmUserProperties
QMap< QString, QList< QVariant > > * _algorithmUserVectorProperties

Additional Inherited Members

Signals inherited from VtkAlgorithmProperties
void ScalarVisibilityChanged (bool on)

Constructor & Destructor Documentation

◆ VtkSurfacesSource()

VtkSurfacesSource::VtkSurfacesSource ( )
protected

Definition at line 25 of file VtkSurfacesSource.cpp.

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}
vtkProperty * GetProperties() const
Returns the vtk properties.
Color getRandomColor()
Returns a random RGB colour.
Definition Color.cpp:18
std::array< unsigned char, 4 > Color
Definition Color.h:12

References VtkAlgorithmProperties::_removable, VtkAlgorithmProperties::GetProperties(), and DataHolderLib::getRandomColor().

Referenced by New(), and vtkTypeMacro().

◆ ~VtkSurfacesSource()

VtkSurfacesSource::~VtkSurfacesSource ( )
overrideprotecteddefault

Member Function Documentation

◆ New()

VtkSurfacesSource * VtkSurfacesSource::New ( )
static

Create new objects with New() because of VTKs object reference counting.

References VtkSurfacesSource().

◆ PrintSelf()

void VtkSurfacesSource::PrintSelf ( ostream & os,
vtkIndent indent )
override

Prints its data on a stream.

Definition at line 37 of file VtkSurfacesSource.cpp.

38{
39 this->Superclass::PrintSelf(os, indent);
40
41 if (_surfaces->empty())
42 {
43 return;
44 }
45
46 os << indent << "== VtkSurfacesSource =="
47 << "\n";
48}
const std::vector< GeoLib::Surface * > * _surfaces
The surfaces to visualize.

References _surfaces.

◆ RequestData()

int VtkSurfacesSource::RequestData ( vtkInformation * request,
vtkInformationVector ** inputVector,
vtkInformationVector * outputVector )
overrideprotected

Computes the polygonal data object.

Definition at line 50 of file VtkSurfacesSource.cpp.

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}

References _surfaces.

◆ RequestInformation()

int VtkSurfacesSource::RequestInformation ( vtkInformation * request,
vtkInformationVector ** inputVector,
vtkInformationVector * outputVector )
overrideprotected

Definition at line 124 of file VtkSurfacesSource.cpp.

128{
129 return 1;
130}

◆ setSurfaces()

void VtkSurfacesSource::setSurfaces ( const std::vector< GeoLib::Surface * > * surfaces)
inline

Sets the surfaces vector.

Definition at line 25 of file VtkSurfacesSource.h.

25{ _surfaces = surfaces; }

References _surfaces.

Referenced by GeoObjectListItem::GeoObjectListItem().

◆ SetUserProperty()

void VtkSurfacesSource::SetUserProperty ( QString name,
QVariant value )
overridevirtual

Generates random colors for each surface.

Reimplemented from VtkAlgorithmProperties.

Definition at line 132 of file VtkSurfacesSource.cpp.

133{
135 (*_algorithmUserProperties)[name] = value;
136}
virtual void SetUserProperty(QString name, QVariant value)
Sets a user property. This should be implemented by subclasses.

References VtkAlgorithmProperties::SetUserProperty().

◆ vtkTypeMacro()

VtkSurfacesSource::vtkTypeMacro ( VtkSurfacesSource ,
vtkPolyDataAlgorithm  )

References VtkSurfacesSource().

Member Data Documentation

◆ _surfaces

const std::vector<GeoLib::Surface*>* VtkSurfacesSource::_surfaces {nullptr}
protected

The surfaces to visualize.

Definition at line 51 of file VtkSurfacesSource.h.

51{nullptr};

Referenced by PrintSelf(), RequestData(), and setSurfaces().


The documentation for this class was generated from the following files: