OGS
VtkStationSource.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// ** VTK INCLUDES **
5#include "VtkStationSource.h"
6
7#include <vtkCellArray.h>
8#include <vtkCellData.h>
9#include <vtkDoubleArray.h>
10#include <vtkInformation.h>
11#include <vtkInformationVector.h>
12#include <vtkLine.h>
13#include <vtkPointData.h>
14#include <vtkPoints.h>
15#include <vtkPolyData.h>
16#include <vtkProperty.h>
17#include <vtkSmartPointer.h>
18#include <vtkStreamingDemandDrivenPipeline.h>
19
20#include "BaseLib/Logging.h"
22#include "vtkObjectFactory.h"
23
25
27{
28 _removable = false; // From VtkAlgorithmProperties
29 this->SetNumberOfInputPorts(0);
30
32 GetProperties()->SetColor(c[0] / 255.0, c[1] / 255.0, c[2] / 255.0);
33}
34
35void VtkStationSource::PrintSelf(ostream& os, vtkIndent indent)
36{
37 this->Superclass::PrintSelf(os, indent);
38
39 if (_stations->empty())
40 {
41 return;
42 }
43
44 os << indent << "== VtkStationSource =="
45 << "\n";
46
47 int i = 0;
48 for (auto station : *_stations)
49 {
50 os << indent << "Station " << i << " (" << (*station)[0] << ", "
51 << (*station)[1] << ", " << (*station)[2] << ")\n";
52 i++;
53 }
54}
55
57int VtkStationSource::RequestData(vtkInformation* request,
58 vtkInformationVector** inputVector,
59 vtkInformationVector* outputVector)
60{
61 (void)request;
62 (void)inputVector;
63
64 if (!_stations)
65 {
66 return 0;
67 }
68 std::size_t nStations = _stations->size();
69 if (nStations == 0)
70 {
71 return 0;
72 }
73
74 bool useStationValues(false);
75 double sValue =
76 static_cast<GeoLib::Station*>((*_stations)[0])->getStationValue();
77 for (std::size_t i = 1; i < nStations; i++)
78 {
79 if (static_cast<GeoLib::Station*>((*_stations)[i])->getStationValue() !=
80 sValue)
81 {
82 useStationValues = true;
83 break;
84 }
85 }
86
87 bool isBorehole = dynamic_cast<GeoLib::StationBorehole*>((*_stations)[0]);
88
89 vtkSmartPointer<vtkInformation> outInfo =
90 outputVector->GetInformationObject(0);
91 vtkSmartPointer<vtkPolyData> output =
92 vtkPolyData::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
93
94 vtkSmartPointer<vtkPoints> newStations = vtkSmartPointer<vtkPoints>::New();
95 vtkSmartPointer<vtkCellArray> newVerts =
96 vtkSmartPointer<vtkCellArray>::New();
97 newVerts->Allocate(nStations);
98
99 vtkSmartPointer<vtkCellArray> newLines;
100
101 if (isBorehole)
102 {
103 newLines = vtkSmartPointer<vtkCellArray>::New();
104 }
105
106 if (outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER()) >
107 0)
108 {
109 return 1;
110 }
111
112 vtkSmartPointer<vtkIntArray> station_ids =
113 vtkSmartPointer<vtkIntArray>::New();
114 station_ids->SetNumberOfComponents(1);
115 station_ids->SetName("SiteIDs");
116
117 vtkSmartPointer<vtkDoubleArray> station_values =
118 vtkSmartPointer<vtkDoubleArray>::New();
119 station_values->SetNumberOfComponents(1);
120 station_values->SetName("StationValue");
121
122 vtkSmartPointer<vtkIntArray> strat_ids =
123 vtkSmartPointer<vtkIntArray>::New();
124 strat_ids->SetNumberOfComponents(1);
125 strat_ids->SetName("Stratigraphies");
126
127 std::size_t lastMaxIndex(0);
128 std::size_t site_count(0);
129
130 // Generate graphic objects
131 for (auto station : *_stations)
132 {
133 vtkIdType sid = newStations->InsertNextPoint(station->data());
134 station_ids->InsertNextValue(site_count);
135 if (useStationValues)
136 {
137 station_values->InsertNextValue(
138 static_cast<GeoLib::Station*>(station)->getStationValue());
139 }
140
141 if (!isBorehole)
142 {
143 newVerts->InsertNextCell(1, &sid);
144 }
145 else
146 {
147 std::vector<GeoLib::Point*> profile =
148 static_cast<GeoLib::StationBorehole*>(station)->getProfile();
149 std::vector<std::string> soilNames =
150 static_cast<GeoLib::StationBorehole*>(station)->getSoilNames();
151 const std::size_t nLayers = profile.size();
152
153 for (std::size_t i = 1; i < nLayers; i++)
154 {
155 newStations->InsertNextPoint(profile[i]->data());
156 station_ids->InsertNextValue(site_count);
157 newLines->InsertNextCell(2);
158 newLines->InsertCellPoint(
159 lastMaxIndex); // start of borehole-layer
160 newLines->InsertCellPoint(
161 ++lastMaxIndex); // end of boreholelayer
162 strat_ids->InsertNextValue(this->GetIndexByName(soilNames[i]));
163 if (useStationValues)
164 {
165 station_values->InsertNextValue(
166 static_cast<GeoLib::Station*>(station)
167 ->getStationValue());
168 }
169 }
170 lastMaxIndex++;
171 }
172 site_count++;
173 }
174
175 output->SetPoints(newStations);
176
177 if (!isBorehole)
178 {
179 output->SetVerts(newVerts);
180 output->GetCellData()->AddArray(station_ids);
181 output->GetCellData()->SetActiveAttribute(
182 "SiteIDs", vtkDataSetAttributes::SCALARS);
183 }
184 else
185 {
186 output->SetLines(newLines);
187 // output->GetCellData()->AddArray(station_ids);
188 output->GetCellData()->AddArray(strat_ids);
189 output->GetCellData()->SetActiveAttribute(
190 "Stratigraphies", vtkDataSetAttributes::SCALARS);
191 }
192 if (useStationValues)
193 {
194 output->GetPointData()->AddArray(station_values);
195 }
196
197 output->Squeeze();
198
199 return 1;
200}
201
202int VtkStationSource::RequestInformation(vtkInformation* /*request*/,
203 vtkInformationVector** /*inputVector*/,
204 vtkInformationVector* /*outputVector*/)
205{
206 return 1;
207}
208
209void VtkStationSource::SetUserProperty(QString name, QVariant value)
210{
211 Q_UNUSED(name);
212 Q_UNUSED(value);
213}
214
215std::size_t VtkStationSource::GetIndexByName(std::string const& name)
216{
217 vtkIdType max_key(0);
218 for (auto& it : _id_map)
219 {
220 if (name == it.first)
221 {
222 return it.second;
223 }
224 if (it.second > max_key)
225 {
226 max_key = it.second;
227 }
228 }
229
230 vtkIdType new_index = (_id_map.empty()) ? 0 : (max_key + 1);
231 INFO("Key '{:s}' has been assigned index {:d}.", name, new_index);
232 _id_map.insert(std::pair<std::string, vtkIdType>(name, new_index));
233 return new_index;
234}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
vtkStandardNewMacro(VtkStationSource)
A borehole as a geometric object.
A Station (observation site) is basically a Point with some additional information.
Definition Station.h:26
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:18
std::array< unsigned char, 4 > Color
Definition Color.h:12