OGS
VtkStationSource Class Reference

Detailed Description

VTK source object for the visualisation of station data (including boreholes)

Definition at line 27 of file VtkStationSource.h.

#include <VtkStationSource.h>

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

Public Member Functions

 vtkTypeMacro (VtkStationSource, vtkPolyDataAlgorithm)
 
const std::map< std::string, DataHolderLib::Color > & getColorLookupTable () const
 
void setStations (const std::vector< GeoLib::Point * > *stations)
 Sets a predefined color lookup table for the colouring of borehole stratigraphies.
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Prints its data on a stream.
 
void SetUserProperty (QString name, QVariant value) override
 Sets a user property. This should be implemented by subclasses.
 
- 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 VtkStationSourceNew ()
 Create new objects with New() because of VTKs object reference counting.
 

Protected Member Functions

 VtkStationSource ()
 
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::Point * > * _stations {nullptr}
 The stations to visualize.
 
std::map< std::string, DataHolderLib::Color_colorLookupTable
 
- 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
 

Private Member Functions

std::size_t GetIndexByName (std::string const &name)
 

Private Attributes

std::map< std::string, vtkIdType > _id_map
 

Additional Inherited Members

- Signals inherited from VtkAlgorithmProperties
void ScalarVisibilityChanged (bool on)
 

Constructor & Destructor Documentation

◆ VtkStationSource()

VtkStationSource::VtkStationSource ( )
protected

Definition at line 37 of file VtkStationSource.cpp.

38{
39 _removable = false; // From VtkAlgorithmProperties
40 this->SetNumberOfInputPorts(0);
41
43 GetProperties()->SetColor(c[0] / 255.0, c[1] / 255.0, c[2] / 255.0);
44}
vtkProperty * GetProperties() const
Returns the vtk properties.
Color getRandomColor()
Returns a random RGB colour.
Definition Color.cpp:29
std::array< unsigned char, 4 > Color
Definition Color.h:24

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

Member Function Documentation

◆ getColorLookupTable()

const std::map< std::string, DataHolderLib::Color > & VtkStationSource::getColorLookupTable ( ) const
inline

Returns the colour lookup table generated for boreholes. This method should only be called after the colour lookup table has actually been build (via RequestData() or setColorLookupTable()).

Definition at line 37 of file VtkStationSource.h.

38 { return _colorLookupTable; }
std::map< std::string, DataHolderLib::Color > _colorLookupTable

References _colorLookupTable.

Referenced by StationTreeView::displayStratigraphy(), and StationTreeView::writeStratigraphiesAsImages().

◆ GetIndexByName()

std::size_t VtkStationSource::GetIndexByName ( std::string const & name)
private

Definition at line 226 of file VtkStationSource.cpp.

227{
228 vtkIdType max_key(0);
229 for (auto& it : _id_map)
230 {
231 if (name == it.first)
232 {
233 return it.second;
234 }
235 if (it.second > max_key)
236 {
237 max_key = it.second;
238 }
239 }
240
241 vtkIdType new_index = (_id_map.empty()) ? 0 : (max_key + 1);
242 INFO("Key '{:s}' has been assigned index {:d}.", name, new_index);
243 _id_map.insert(std::pair<std::string, vtkIdType>(name, new_index));
244 return new_index;
245}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
std::map< std::string, vtkIdType > _id_map

References _id_map, and INFO().

Referenced by RequestData().

◆ New()

static VtkStationSource * VtkStationSource::New ( )
static

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

◆ PrintSelf()

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

Prints its data on a stream.

Definition at line 46 of file VtkStationSource.cpp.

47{
48 this->Superclass::PrintSelf(os, indent);
49
50 if (_stations->empty())
51 {
52 return;
53 }
54
55 os << indent << "== VtkStationSource =="
56 << "\n";
57
58 int i = 0;
59 for (auto station : *_stations)
60 {
61 os << indent << "Station " << i << " (" << (*station)[0] << ", "
62 << (*station)[1] << ", " << (*station)[2] << ")\n";
63 i++;
64 }
65}
const std::vector< GeoLib::Point * > * _stations
The stations to visualize.

References _stations.

◆ RequestData()

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

Computes the polygonal data object.

Create 3d Station objects.

Definition at line 68 of file VtkStationSource.cpp.

71{
72 (void)request;
73 (void)inputVector;
74
75 if (!_stations)
76 {
77 return 0;
78 }
79 std::size_t nStations = _stations->size();
80 if (nStations == 0)
81 {
82 return 0;
83 }
84
85 bool useStationValues(false);
86 double sValue =
87 static_cast<GeoLib::Station*>((*_stations)[0])->getStationValue();
88 for (std::size_t i = 1; i < nStations; i++)
89 {
90 if (static_cast<GeoLib::Station*>((*_stations)[i])->getStationValue() !=
91 sValue)
92 {
93 useStationValues = true;
94 break;
95 }
96 }
97
98 bool isBorehole = dynamic_cast<GeoLib::StationBorehole*>((*_stations)[0]);
99
100 vtkSmartPointer<vtkInformation> outInfo =
101 outputVector->GetInformationObject(0);
102 vtkSmartPointer<vtkPolyData> output =
103 vtkPolyData::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
104
105 vtkSmartPointer<vtkPoints> newStations = vtkSmartPointer<vtkPoints>::New();
106 vtkSmartPointer<vtkCellArray> newVerts =
107 vtkSmartPointer<vtkCellArray>::New();
108 newVerts->Allocate(nStations);
109
110 vtkSmartPointer<vtkCellArray> newLines;
111
112 if (isBorehole)
113 {
114 newLines = vtkSmartPointer<vtkCellArray>::New();
115 }
116
117 if (outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER()) >
118 0)
119 {
120 return 1;
121 }
122
123 vtkSmartPointer<vtkIntArray> station_ids =
124 vtkSmartPointer<vtkIntArray>::New();
125 station_ids->SetNumberOfComponents(1);
126 station_ids->SetName("SiteIDs");
127
128 vtkSmartPointer<vtkDoubleArray> station_values =
129 vtkSmartPointer<vtkDoubleArray>::New();
130 station_values->SetNumberOfComponents(1);
131 station_values->SetName("StationValue");
132
133 vtkSmartPointer<vtkIntArray> strat_ids =
134 vtkSmartPointer<vtkIntArray>::New();
135 strat_ids->SetNumberOfComponents(1);
136 strat_ids->SetName("Stratigraphies");
137
138 std::size_t lastMaxIndex(0);
139 std::size_t site_count(0);
140
141 // Generate graphic objects
142 for (auto station : *_stations)
143 {
144 vtkIdType sid = newStations->InsertNextPoint(station->data());
145 station_ids->InsertNextValue(site_count);
146 if (useStationValues)
147 {
148 station_values->InsertNextValue(
149 static_cast<GeoLib::Station*>(station)->getStationValue());
150 }
151
152 if (!isBorehole)
153 {
154 newVerts->InsertNextCell(1, &sid);
155 }
156 else
157 {
158 std::vector<GeoLib::Point*> profile =
159 static_cast<GeoLib::StationBorehole*>(station)->getProfile();
160 std::vector<std::string> soilNames =
161 static_cast<GeoLib::StationBorehole*>(station)->getSoilNames();
162 const std::size_t nLayers = profile.size();
163
164 for (std::size_t i = 1; i < nLayers; i++)
165 {
166 newStations->InsertNextPoint(profile[i]->data());
167 station_ids->InsertNextValue(site_count);
168 newLines->InsertNextCell(2);
169 newLines->InsertCellPoint(
170 lastMaxIndex); // start of borehole-layer
171 newLines->InsertCellPoint(
172 ++lastMaxIndex); // end of boreholelayer
173 strat_ids->InsertNextValue(this->GetIndexByName(soilNames[i]));
174 if (useStationValues)
175 {
176 station_values->InsertNextValue(
177 static_cast<GeoLib::Station*>(station)
178 ->getStationValue());
179 }
180 }
181 lastMaxIndex++;
182 }
183 site_count++;
184 }
185
186 output->SetPoints(newStations);
187
188 if (!isBorehole)
189 {
190 output->SetVerts(newVerts);
191 output->GetCellData()->AddArray(station_ids);
192 output->GetCellData()->SetActiveAttribute(
193 "SiteIDs", vtkDataSetAttributes::SCALARS);
194 }
195 else
196 {
197 output->SetLines(newLines);
198 // output->GetCellData()->AddArray(station_ids);
199 output->GetCellData()->AddArray(strat_ids);
200 output->GetCellData()->SetActiveAttribute(
201 "Stratigraphies", vtkDataSetAttributes::SCALARS);
202 }
203 if (useStationValues)
204 {
205 output->GetPointData()->AddArray(station_values);
206 }
207
208 output->Squeeze();
209
210 return 1;
211}
A borehole as a geometric object.
A Station (observation site) is basically a Point with some additional information.
Definition Station.h:37
double getStationValue() const
Returns the specific value for this station.
Definition Station.h:73
std::size_t GetIndexByName(std::string const &name)
bool isBorehole(GeoLib::Point const *pnt)

References _stations, GetIndexByName(), and GeoLib::Station::getStationValue().

◆ RequestInformation()

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

Definition at line 213 of file VtkStationSource.cpp.

216{
217 return 1;
218}

◆ setStations()

void VtkStationSource::setStations ( const std::vector< GeoLib::Point * > * stations)
inline

Sets a predefined color lookup table for the colouring of borehole stratigraphies.

Sets the stations as a vector

Definition at line 45 of file VtkStationSource.h.

45{ _stations = stations; }

References _stations.

◆ SetUserProperty()

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

Sets a user property. This should be implemented by subclasses.

Reimplemented from VtkAlgorithmProperties.

Definition at line 220 of file VtkStationSource.cpp.

221{
222 Q_UNUSED(name);
223 Q_UNUSED(value);
224}

◆ vtkTypeMacro()

VtkStationSource::vtkTypeMacro ( VtkStationSource ,
vtkPolyDataAlgorithm  )

Member Data Documentation

◆ _colorLookupTable

std::map<std::string, DataHolderLib::Color> VtkStationSource::_colorLookupTable
protected

The colour table for stratigraphic data. This table is either set using the setColorLookupTable() method or is generated automatically with random colours while creating the VtkStationSource-object.

Definition at line 69 of file VtkStationSource.h.

Referenced by getColorLookupTable().

◆ _id_map

std::map<std::string, vtkIdType> VtkStationSource::_id_map
private

Definition at line 74 of file VtkStationSource.h.

Referenced by GetIndexByName().

◆ _stations

const std::vector<GeoLib::Point*>* VtkStationSource::_stations {nullptr}
protected

The stations to visualize.

Definition at line 65 of file VtkStationSource.h.

65{nullptr};

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


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