OGS
VtkGeoImageSource Class Reference

Detailed Description

The VtkVisPipeline source object of a geo-referenced image (file).

Definition at line 39 of file VtkGeoImageSource.h.

#include <VtkGeoImageSource.h>

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

Public Member Functions

 vtkTypeMacro (VtkGeoImageSource, vtkSimpleImageToImageFilter)
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Prints information about itself.
 
 VtkGeoImageSource (const VtkGeoImageSource &)=delete
 
void operator= (const VtkGeoImageSource &)=delete
 
vtkImageData * getImageData ()
 Returns the ImageData object.
 
bool readImage (const QString &filename)
 Reads an image from file.
 
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.
 
- 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 VtkGeoImageSourceNew ()
 Create new objects with New() because of VTKs reference counting.
 
static std::optional< GeoLib::RasterconvertToRaster (VtkGeoImageSource *const source)
 

Protected Member Functions

 VtkGeoImageSource ()
 Constructor.
 
 ~VtkGeoImageSource () override
 Destructor.
 
void SimpleExecute (vtkImageData *input, vtkImageData *output) override
 Filter execution.
 

Private Attributes

vtkImageAlgorithm * _imageSource {nullptr}
 

Additional Inherited Members

- Signals inherited from VtkAlgorithmProperties
void ScalarVisibilityChanged (bool on)
 
- 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
 

Constructor & Destructor Documentation

◆ VtkGeoImageSource() [1/2]

VtkGeoImageSource::VtkGeoImageSource ( const VtkGeoImageSource & )
delete

◆ VtkGeoImageSource() [2/2]

VtkGeoImageSource::VtkGeoImageSource ( )
protecteddefault

Constructor.

◆ ~VtkGeoImageSource()

VtkGeoImageSource::~VtkGeoImageSource ( )
overrideprotected

Destructor.

Definition at line 63 of file VtkGeoImageSource.cpp.

64{
65 if (_imageSource)
66 {
67 _imageSource->Delete();
68 }
69}
vtkImageAlgorithm * _imageSource

References _imageSource.

Member Function Documentation

◆ convertToRaster()

std::optional< GeoLib::Raster > VtkGeoImageSource::convertToRaster ( VtkGeoImageSource *const source)
static

Definition at line 122 of file VtkGeoImageSource.cpp.

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
Contains the relevant information when storing a geoscientific raster data.
Definition Raster.h:28

References ERR().

Referenced by VtkVisPipelineView::showImageToMeshConversionDialog(), and VtkVisImageItem::writeAsRaster().

◆ getImageData()

vtkImageData * VtkGeoImageSource::getImageData ( )

Returns the ImageData object.

Definition at line 94 of file VtkGeoImageSource.cpp.

95{
96 return this->_imageSource->GetImageDataInput(0);
97}

References _imageSource.

◆ New()

static VtkGeoImageSource * VtkGeoImageSource::New ( )
static

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

Referenced by NetCdfConfigureDialog::createDataObject(), and MainWindow::loadFile().

◆ operator=()

void VtkGeoImageSource::operator= ( const VtkGeoImageSource & )
delete

◆ PrintSelf()

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

Prints information about itself.

Definition at line 71 of file VtkGeoImageSource.cpp.

72{
73 this->Superclass::PrintSelf(os, indent);
74}

◆ readImage()

bool VtkGeoImageSource::readImage ( const QString & filename)

Reads an image from file.

Definition at line 76 of file VtkGeoImageSource.cpp.

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}
void setImage(vtkImageAlgorithm *image, const QString &name)
Imports an existing image object.
static vtkImageAlgorithm * loadImage(const std::string &fileName)
Loads an image- or raster-file into an vtkImageAlgorithm-Object.
Definition VtkRaster.cpp:43

References VtkRaster::loadImage(), and setImage().

Referenced by MainWindow::loadFile().

◆ setImage()

void VtkGeoImageSource::setImage ( vtkImageAlgorithm * image,
const QString & name )

Imports an existing image object.

Definition at line 87 of file VtkGeoImageSource.cpp.

88{
89 this->_imageSource = image;
90 this->SetInputConnection(_imageSource->GetOutputPort());
91 this->SetName(name);
92}
void SetName(QString name)
Sets the name.

References _imageSource, and VtkAlgorithmProperties::SetName().

Referenced by NetCdfConfigureDialog::createDataObject(), and readImage().

◆ SetUserProperty()

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

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

Reimplemented from VtkAlgorithmProperties.

Definition at line 116 of file VtkGeoImageSource.cpp.

117{
118 Q_UNUSED(name);
119 Q_UNUSED(value);
120}

◆ SimpleExecute()

void VtkGeoImageSource::SimpleExecute ( vtkImageData * input,
vtkImageData * output )
overrideprotected

Filter execution.

Definition at line 99 of file VtkGeoImageSource.cpp.

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}
void vtkSimpleImageFilterExampleExecute(vtkImageData *input, vtkImageData *output, IT *inPtr, IT *outPtr)

References vtkSimpleImageFilterExampleExecute().

◆ vtkTypeMacro()

VtkGeoImageSource::vtkTypeMacro ( VtkGeoImageSource ,
vtkSimpleImageToImageFilter  )

Member Data Documentation

◆ _imageSource

vtkImageAlgorithm* VtkGeoImageSource::_imageSource {nullptr}
private

Definition at line 79 of file VtkGeoImageSource.h.

79{nullptr};

Referenced by ~VtkGeoImageSource(), getImageData(), and setImage().


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