OGS
VtkVisImageItem.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// ** INCLUDES **
5#include "VtkVisImageItem.h"
6
7#include <vtkActor.h>
8#include <vtkDataArray.h>
9#include <vtkDataSetMapper.h>
10#include <vtkImageAlgorithm.h>
11#include <vtkImageChangeInformation.h>
12#include <vtkImageData.h>
13#include <vtkImageShiftScale.h>
14#include <vtkPointData.h>
15#include <vtkRenderer.h>
16#include <vtkSmartPointer.h>
17
18#include <QFileDialog>
19
21#include "Base/OGSError.h"
22#include "BaseLib/FileTools.h"
23#include "BaseLib/Logging.h"
25#include "GeoLib/Raster.h"
27#include "VtkGeoImageSource.h"
28
29// export
30#include <vtkImageActor.h>
31#include <vtkXMLImageDataWriter.h>
32
34 vtkAlgorithm* algorithm, TreeItem* parentItem,
35 const QList<QVariant> data /*= QList<QVariant>()*/)
37{
38}
39
47
52
53void VtkVisImageItem::Initialize(vtkRenderer* renderer)
54{
55 auto* image = dynamic_cast<vtkImageAlgorithm*>(_algorithm);
56 if (!image)
57 {
58 OGSError::box("Cast to image failed.");
59 return;
60 }
61 image->Update();
62
63 double origin[3];
64 double spacing[3];
65 double range[2];
66 vtkImageData* image_data = image->GetOutput();
67 if (!image_data)
68 {
69 OGSError::box("GetOutput() - Image could not be initialized.");
70 return;
71 }
72 image_data->GetOrigin(origin);
73 image_data->GetSpacing(spacing);
74 vtkPointData* point_data = image_data->GetPointData();
75 if (point_data)
76 {
77 auto* scalars = point_data->GetScalars();
78 if (scalars)
79 {
80 scalars->GetRange(range);
81 }
82 }
83 vtkImageShiftScale* scale = vtkImageShiftScale::New();
84 scale->SetOutputScalarTypeToUnsignedChar();
85 scale->SetInputConnection(image->GetOutputPort());
86 scale->SetShift(-range[0]);
87 scale->SetScale(255.0 / (range[1] - range[0]));
88
89 _transformFilter = vtkImageChangeInformation::New();
90 _transformFilter->SetInputConnection(scale->GetOutputPort());
91 // double origin[3];
92 // img->getOrigin(origin);
93 // double spacing = img->getSpacing();
94 //_transformFilter->SetOutputOrigin(origin2);
95 //_transformFilter->SetOutputSpacing(spacing2);
96 _transformFilter->Update();
97
98 _renderer = renderer;
99
100 // Use a special vtkImageActor instead of vtkActor
101 vtkImageActor* imageActor = vtkImageActor::New();
102 imageActor->SetInputData(_transformFilter->GetOutput());
103 _actor = imageActor;
104 _renderer->AddActor(_actor);
105
106 // Set pre-set properties
107 auto* vtkProps = dynamic_cast<VtkAlgorithmProperties*>(_algorithm);
108 if (vtkProps)
109 {
110 setVtkProperties(vtkProps);
111 }
112
113 auto* parentItem = dynamic_cast<VtkVisPipelineItem*>(this->parentItem());
114 while (parentItem)
115 {
116 auto* parentProps =
117 dynamic_cast<VtkAlgorithmProperties*>(parentItem->algorithm());
118 if (parentProps)
119 {
120 auto* newProps = new VtkAlgorithmProperties();
121 newProps->SetScalarVisibility(parentProps->GetScalarVisibility());
122 newProps->SetTexture(parentProps->GetTexture());
123 setVtkProperties(newProps);
124 vtkProps = newProps;
125 parentItem = nullptr;
126 }
127 else
128 {
129 parentItem =
130 dynamic_cast<VtkVisPipelineItem*>(parentItem->parentItem());
131 }
132 }
133
134 // Set active scalar to the desired one from VtkAlgorithmProperties
135 // or to match those of the parent.
136 if (vtkProps)
137 {
138 if (vtkProps->GetActiveAttribute().length() > 0)
139 {
140 this->SetActiveAttribute(vtkProps->GetActiveAttribute());
141 }
142 else
143 {
144 auto* visParentItem =
145 dynamic_cast<VtkVisPipelineItem*>(this->parentItem());
146 if (visParentItem)
147 {
148 this->SetActiveAttribute(visParentItem->GetActiveAttribute());
149 }
150 if (vtkProps->GetTexture() != nullptr)
151 {
152 this->SetActiveAttribute("Solid Color");
153 }
154 }
155 }
156}
157
159{
160 // todo implementation
161 (void)vtkProps;
162}
163
165 const std::string& filename) const
166{
167 std::string file_name_cpy(filename);
168 auto* algID = dynamic_cast<vtkImageAlgorithm*>(algorithm);
169 if (algID)
170 {
171 vtkSmartPointer<vtkXMLImageDataWriter> iWriter =
172 vtkSmartPointer<vtkXMLImageDataWriter>::New();
173 iWriter->SetInputData(algID->GetOutputDataObject(0));
174 if (BaseLib::getFileExtension(filename) != ".vti")
175 {
176 file_name_cpy.append(".vti");
177 }
178 iWriter->SetFileName(file_name_cpy.c_str());
179 return iWriter->Write();
180 }
181 ERR("VtkVisPipelineItem::writeToFile() - Unknown data type.");
182 return 0;
183}
184
185void VtkVisImageItem::setTranslation(double x, double y, double z) const
186{
187 _transformFilter->SetOriginTranslation(x, y, z);
188}
189
191{
192 return _transformFilter;
193}
194
196{
197 vtkSmartPointer<VtkGeoImageSource> imageSource =
198 VtkGeoImageSource::SafeDownCast(_algorithm);
199
200 auto const raster = VtkGeoImageSource::convertToRaster(imageSource);
201
202 if (!raster)
203 {
204 OGSError::box("Image could not be converted to a raster.");
205 return false;
206 }
207
208 QFileInfo const info(this->data(0).toString());
209 QString const rastername = info.completeBaseName() + ".asc";
210 QString const filetype("ESRI ASCII raster file (*.asc)");
211 QString const filename = QFileDialog::getSaveFileName(
212 nullptr, "Save raster as",
213 LastSavedFileDirectory::getDir() + rastername, filetype);
216 filename.toStdString());
217 return true;
218}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
static void writeRasterAsASC(GeoLib::Raster const &raster, std::string const &file_name)
Writes an Esri asc-file.
static void setDir(const QString &path)
Sets the directory last used for saving a file.
static const QString getDir()
Returns the directory last used for saving a file.
static void box(const QString &e)
Definition OGSError.cpp:13
TreeItem * parentItem() const
Definition TreeItem.cpp:104
TreeItem(QList< QVariant > data, TreeItem *parent)
Definition TreeItem.cpp:12
Contains properties for the visualization of objects as VtkVisPipelineItems.
Is used to combine several filter in one VtkVisPipelineItem. You can use vtk filter and custom filter...
static std::optional< GeoLib::Raster > convertToRaster(VtkGeoImageSource *const source)
vtkImageChangeInformation * _transformFilter
VtkVisImageItem(vtkAlgorithm *algorithm, TreeItem *parentItem, const QList< QVariant > data=QList< QVariant >())
Constructor for a source/filter object.
bool writeAsRaster()
Allows writing this item's source object as an ASCII raster file.
void Initialize(vtkRenderer *renderer) override
Initializes vtkMapper and vtkActor necessary for visualization of the item and sets the item's proper...
vtkAlgorithm * transformFilter() const override
int callVTKWriter(vtkAlgorithm *algorithm, const std::string &filename) const override
Selects the appropriate VTK-Writer object and writes the object to a file with the given name.
~VtkVisImageItem() override
void setVtkProperties(VtkAlgorithmProperties *vtkProps)
void setTranslation(double x, double y, double z) const override
Translates the item in visualisation-space. This function is empty and needs to be implemented by der...
vtkAlgorithm * algorithm() const
Returns the algorithm object.
VtkVisPipelineItem(vtkAlgorithm *algorithm, TreeItem *parentItem, const QList< QVariant > data=QList< QVariant >())
Constructor for a source/filter object.
vtkAlgorithm * _algorithm
QVariant data(int column) const override
VtkCompositeFilter * compositeFilter() const
Returns the composite filter.
virtual void SetActiveAttribute(const QString &str)
std::string getFileExtension(const std::string &path)