OGS
VtkGeoImageSource.cpp
Go to the documentation of this file.
1
15// ** INCLUDES **
16#include "VtkGeoImageSource.h"
17
18#include <vtkFloatArray.h>
19#include <vtkImageChangeInformation.h>
20#include <vtkImageData.h>
21#include <vtkImageImport.h>
22#include <vtkImageShiftScale.h>
23#include <vtkIntArray.h>
24#include <vtkObjectFactory.h>
25#include <vtkPointData.h>
26
27#include "VtkRaster.h"
28
30
31// This function is copied from the vtkSimpleImageFilterExample.cpp
32//
33// The switch statement in Execute will call this method with
34// the appropriate input type (IT). Note that this example assumes
35// that the output data type is the same as the input data type.
36// This is not always the case.
37template <class IT>
38void vtkSimpleImageFilterExampleExecute(vtkImageData* input,
39 vtkImageData* output, IT* inPtr,
40 IT* outPtr)
41{
42 int dims[3];
43 input->GetDimensions(dims);
44 if (input->GetScalarType() != output->GetScalarType())
45 {
46 vtkGenericWarningMacro(
47 << "Execute: input ScalarType, " << input->GetScalarType()
48 << ", must match out ScalarType " << output->GetScalarType());
49 return;
50 }
51 // HACK LB Multiply by number of scalar components due to RGBA values ?????
52 int size =
53 dims[0] * dims[1] * dims[2] * input->GetNumberOfScalarComponents();
54
55 for (int i = 0; i < size; i++)
56 {
57 outPtr[i] = inPtr[i];
58 }
59}
60
62
64{
65 if (_imageSource)
66 {
67 _imageSource->Delete();
68 }
69}
70
71void VtkGeoImageSource::PrintSelf(ostream& os, vtkIndent indent)
72{
73 this->Superclass::PrintSelf(os, indent);
74}
75
76bool VtkGeoImageSource::readImage(const QString& filename)
77{
78 vtkImageAlgorithm* img(VtkRaster::loadImage(filename.toStdString()));
79 if (img == nullptr)
80 {
81 return false;
82 }
83 this->setImage(img, filename);
84 return true;
85}
86
87void VtkGeoImageSource::setImage(vtkImageAlgorithm* image, const QString& name)
88{
89 this->_imageSource = image;
90 this->SetInputConnection(_imageSource->GetOutputPort());
91 this->SetName(name);
92}
93
95{
96 return this->_imageSource->GetImageDataInput(0);
97}
98
99void VtkGeoImageSource::SimpleExecute(vtkImageData* input, vtkImageData* output)
100{
101 vtkDebugMacro(<< "Executing VtkGeoImageSource");
102 void* inPtr = input->GetScalarPointer();
103 void* outPtr = output->GetScalarPointer();
104 switch (output->GetScalarType())
105 {
106 // This is simply a #define for a big case list.
107 // It handles all data types that VTK supports.
108 vtkTemplateMacro(vtkSimpleImageFilterExampleExecute(
109 input, output, (VTK_TT*)(inPtr), (VTK_TT*)(outPtr)));
110 default:
111 vtkGenericWarningMacro("Execute: Unknown input ScalarType");
112 return;
113 }
114}
115
116void VtkGeoImageSource::SetUserProperty(QString name, QVariant value)
117{
118 Q_UNUSED(name);
119 Q_UNUSED(value);
120}
121
122std::optional<GeoLib::Raster> VtkGeoImageSource::convertToRaster(
123 VtkGeoImageSource* const source)
124{
125 int dims[3];
126 source->GetOutput()->GetDimensions(dims);
127 double origin[3];
128 source->GetOutput()->GetOrigin(origin);
129 double spacing[3];
130 source->GetOutput()->GetSpacing(spacing);
131 MathLib::Point3d const origin_pnt(
132 std::array<double, 3>{{origin[0] - 0.5 * spacing[0],
133 origin[1] - 0.5 * spacing[1], origin[2]}});
134 GeoLib::RasterHeader const header = {static_cast<std::size_t>(dims[0]),
135 static_cast<std::size_t>(dims[1]),
136 static_cast<std::size_t>(dims[2]),
137 origin_pnt,
138 spacing[0],
139 -9999};
140
141 vtkSmartPointer<vtkDataArray> const pixelData =
142 vtkSmartPointer<vtkDataArray>(
143 source->GetOutput()->GetPointData()->GetScalars());
144 int const nTuple = pixelData->GetNumberOfComponents();
145 if (nTuple < 1 || nTuple > 4)
146 {
147 ERR("VtkMeshConverter::convertImgToMesh(): Unsupported pixel "
148 "composition!");
149 return {};
150 }
151
152 std::vector<double> pix(header.n_cols * header.n_rows * header.n_depth, 0);
153 for (std::size_t k = 0; k < header.n_depth; k++)
154 {
155 std::size_t const layer_idx = (k * header.n_rows * header.n_cols);
156 for (std::size_t i = 0; i < header.n_rows; i++)
157 {
158 std::size_t const idx = i * header.n_cols + layer_idx;
159 for (std::size_t j = 0; j < header.n_cols; j++)
160 {
161 double const* const colour = pixelData->GetTuple(idx + j);
162 bool const visible = (nTuple == 2 || nTuple == 4)
163 ? (colour[nTuple - 1] != 0)
164 : true;
165 if (!visible)
166 {
167 pix[idx + j] = header.no_data;
168 }
169 else
170 {
171 pix[idx + j] = (nTuple < 3) ? colour[0] : // grey (+ alpha)
172 (0.3 * colour[0] + 0.6 * colour[1] +
173 0.1 * colour[2]); // rgb(a)
174 }
175 }
176 }
177 }
178
179 return std::make_optional<GeoLib::Raster>(header, pix.begin(), pix.end());
180}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
vtkStandardNewMacro(VtkGeoImageSource)
void vtkSimpleImageFilterExampleExecute(vtkImageData *input, vtkImageData *output, IT *inPtr, IT *outPtr)
Definition of the VtkGeoImageSource class.
Definition of the VtkRaster class.
void SetName(QString name)
Sets the name.
The VtkVisPipeline source object of a geo-referenced image (file).
static std::optional< GeoLib::Raster > convertToRaster(VtkGeoImageSource *const source)
VtkGeoImageSource()
Constructor.
void setImage(vtkImageAlgorithm *image, const QString &name)
Imports an existing image object.
void SetUserProperty(QString name, QVariant value) override
Sets a user property. This should be implemented by subclasses.
vtkImageAlgorithm * _imageSource
~VtkGeoImageSource() override
Destructor.
vtkImageData * getImageData()
Returns the ImageData object.
void PrintSelf(ostream &os, vtkIndent indent) override
Prints information about itself.
void SimpleExecute(vtkImageData *input, vtkImageData *output) override
Filter execution.
bool readImage(const QString &filename)
Reads an image from file.
static vtkImageAlgorithm * loadImage(const std::string &fileName)
Loads an image- or raster-file into an vtkImageAlgorithm-Object.
Definition VtkRaster.cpp:43
Contains the relevant information when storing a geoscientific raster data.
Definition Raster.h:28