OGS
VtkStationSource.cpp
Go to the documentation of this file.
1
15// ** VTK INCLUDES **
16#include "VtkStationSource.h"
17
18#include <vtkCellArray.h>
19#include <vtkCellData.h>
20#include <vtkDoubleArray.h>
21#include <vtkInformation.h>
22#include <vtkInformationVector.h>
23#include <vtkLine.h>
24#include <vtkPointData.h>
25#include <vtkPoints.h>
26#include <vtkPolyData.h>
27#include <vtkProperty.h>
28#include <vtkSmartPointer.h>
29#include <vtkStreamingDemandDrivenPipeline.h>
30
31#include "BaseLib/Logging.h"
33#include "vtkObjectFactory.h"
34
36
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}
45
46void VtkStationSource::PrintSelf(ostream& os, vtkIndent indent)
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}
66
68int VtkStationSource::RequestData(vtkInformation* request,
69 vtkInformationVector** inputVector,
70 vtkInformationVector* outputVector)
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}
212
213int VtkStationSource::RequestInformation(vtkInformation* /*request*/,
214 vtkInformationVector** /*inputVector*/,
215 vtkInformationVector* /*outputVector*/)
216{
217 return 1;
218}
219
220void VtkStationSource::SetUserProperty(QString name, QVariant value)
221{
222 Q_UNUSED(name);
223 Q_UNUSED(value);
224}
225
226std::size_t VtkStationSource::GetIndexByName(std::string const& name)
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
Definition of the StationBorehole class.
vtkStandardNewMacro(VtkStationSource)
Definition of the VtkStationSource class.
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
vtkProperty * GetProperties() const
Returns the vtk properties.
VTK source object for the visualisation of station data (including boreholes)
std::size_t GetIndexByName(std::string const &name)
void SetUserProperty(QString name, QVariant value) override
Sets a user property. This should be implemented by subclasses.
std::map< std::string, vtkIdType > _id_map
const std::vector< GeoLib::Point * > * _stations
The stations to visualize.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Computes the polygonal data object.
int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
void PrintSelf(ostream &os, vtkIndent indent) override
Prints its data on a stream.
Color getRandomColor()
Returns a random RGB colour.
Definition Color.cpp:29
std::array< unsigned char, 4 > Color
Definition Color.h:24