OGS
VtkVisPipelineView.cpp
Go to the documentation of this file.
1
15// ** INCLUDES **
16#include "VtkVisPipelineView.h"
17
18#include <vtkDataSetMapper.h>
19#include <vtkProp3D.h>
20
21#include <QAbstractItemModel>
22#include <QContextMenuEvent>
23#include <QFileDialog>
24#include <QHeaderView>
25#include <QMenu>
26#include <QMessageBox>
27#include <QSettings>
28
30#include "Base/OGSError.h"
32#include "MeshLib/Mesh.h"
34#include "VtkVisImageItem.h"
35#include "VtkVisPipeline.h"
36#include "VtkVisPipelineItem.h"
37#include "VtkVisPointSetItem.h"
38
39// image to mesh conversion
40#include <vtkDataObject.h>
41#include <vtkGenericDataObjectReader.h>
42#include <vtkImageData.h>
43#include <vtkSmartPointer.h>
44#include <vtkTransformFilter.h>
45#include <vtkUnstructuredGrid.h>
46#include <vtkUnstructuredGridAlgorithm.h>
47#include <vtkXMLUnstructuredGridReader.h>
48
50#include "VtkGeoImageSource.h"
51
53 : QTreeView(parent)
54{
55 this->setItemsExpandable(false);
56 auto* checkboxDelegate = new CheckboxDelegate(this);
57 this->setItemDelegateForColumn(1, checkboxDelegate);
58 this->header()->setStretchLastSection(true);
59 this->header()->setSectionResizeMode(QHeaderView::ResizeToContents);
60}
61
62void VtkVisPipelineView::setModel(QAbstractItemModel* model)
63{
64 QTreeView::setModel(model);
65
66 // Move Visible checkbox to the left.
67 // This is done here because at constructor time there aren't any sections.
68 this->header()->moveSection(1, 0);
69}
70
71void VtkVisPipelineView::contextMenuEvent(QContextMenuEvent* event)
72{
73 QModelIndex index = selectionModel()->currentIndex();
74 if (index.isValid())
75 {
76 // check object type
77 VtkVisPipelineItem* item = static_cast<VtkVisPipelineItem*>(
78 static_cast<VtkVisPipeline*>(this->model())
79 ->getItem(this->selectionModel()->currentIndex()));
80 int objectType =
81 item->algorithm()->GetOutputDataObject(0)->GetDataObjectType();
82 VtkAlgorithmProperties* vtkProps = item->getVtkProperties();
83 bool isSourceItem =
84 !(this->selectionModel()->currentIndex().parent().isValid());
85
86 QMenu menu;
87 QAction* addFilterAction = menu.addAction("Add filter...");
88
89 if (objectType == VTK_IMAGE_DATA)
90 {
91 // this exception is needed as image object are only displayed in
92 // the vis-pipeline
93 isSourceItem = false;
94 QAction* saveRasterAction = menu.addAction("Save Raster...");
95 connect(saveRasterAction, SIGNAL(triggered()), this,
96 SLOT(writeRaster()));
97 QAction* addMeshingAction =
98 menu.addAction("Convert Raster to Mesh...");
99 connect(addMeshingAction, SIGNAL(triggered()), this,
101 }
102 else
103 {
104 QAction* addLUTAction = menu.addAction("Add color table...");
105 connect(addLUTAction, SIGNAL(triggered()), this,
106 SLOT(addColorTable()));
107 }
108
109 if (objectType == VTK_UNSTRUCTURED_GRID)
110 {
111 QAction* addConvertToMeshAction =
112 menu.addAction("Convert to Mesh...");
113 connect(addConvertToMeshAction, SIGNAL(triggered()), this,
114 SLOT(convertVTKToOGSMesh()));
115 }
116 menu.addSeparator();
117 QAction* exportVtkAction = menu.addAction("Export as VTK");
118
119 if (!isSourceItem || vtkProps->IsRemovable())
120 {
121 menu.addSeparator();
122 QAction* removeAction = menu.addAction("Remove");
123 connect(removeAction, SIGNAL(triggered()), this,
125 }
126
127 connect(addFilterAction, SIGNAL(triggered()), this,
128 SLOT(addPipelineFilterItem()));
129 connect(exportVtkAction, SIGNAL(triggered()), this,
131
132 menu.exec(event->globalPos());
133 }
134}
135
137{
138 QSettings settings;
139 QModelIndex idx = this->selectionModel()->currentIndex();
140 QString filename = QFileDialog::getSaveFileName(
141 this, "Export object to vtk-file",
142 settings.value("lastExportedFileDirectory").toString(),
143 "VTK file (*.*)");
144 if (!filename.isEmpty())
145 {
146 static_cast<VtkVisPipelineItem*>(
147 static_cast<VtkVisPipeline*>(this->model())->getItem(idx))
148 ->writeToFile(filename.toStdString());
149 QDir dir = QDir(filename);
150 settings.setValue("lastExportedFileDirectory", dir.absolutePath());
151 }
152}
153
155{
156 emit requestRemovePipelineItem(selectionModel()->currentIndex());
157}
158
160{
161 emit requestAddPipelineFilterItem(selectionModel()->currentIndex());
162}
163
165{
167 if (dlg.exec() != QDialog::Accepted)
168 {
169 return;
170 }
171
172 vtkSmartPointer<vtkAlgorithm> algorithm =
173 static_cast<VtkVisPipelineItem*>(
174 static_cast<VtkVisPipeline*>(this->model())
175 ->getItem(this->selectionModel()->currentIndex()))
176 ->algorithm();
177
178 vtkSmartPointer<VtkGeoImageSource> imageSource =
179 VtkGeoImageSource::SafeDownCast(algorithm);
180
181 auto const raster = VtkGeoImageSource::convertToRaster(imageSource);
182 if (!raster)
183 {
184 OGSError::box("Image could not be converted to a raster.");
185 return;
186 }
187
189 *raster, dlg.getElementSelection(), dlg.getIntensitySelection(),
190 dlg.getArrayName());
191 if (mesh)
192 {
193 mesh->setName(dlg.getMeshName());
194 emit meshAdded(mesh.release());
195 }
196 else
197 {
198 OGSError::box("Error creating mesh.");
199 }
200}
201
203{
204 QModelIndex const index = this->selectionModel()->currentIndex();
205 static_cast<VtkVisImageItem*>(
206 static_cast<VtkVisPipeline*>(this->model())->getItem(index))
207 ->writeAsRaster();
208}
209
211{
212 VtkVisPipelineItem* item = static_cast<VtkVisPipelineItem*>(
213 static_cast<VtkVisPipeline*>(this->model())
214 ->getItem(this->selectionModel()->currentIndex()));
215 vtkSmartPointer<vtkAlgorithm> algorithm = item->algorithm();
216
217 vtkUnstructuredGrid* grid(nullptr);
218 vtkUnstructuredGridAlgorithm* ugAlg =
219 vtkUnstructuredGridAlgorithm::SafeDownCast(algorithm);
220 if (ugAlg)
221 {
222 grid = ugAlg->GetOutput();
223 }
224 else
225 {
226 // for old filetypes
227 vtkGenericDataObjectReader* dataReader =
228 vtkGenericDataObjectReader::SafeDownCast(algorithm);
229 if (dataReader)
230 {
231 grid = vtkUnstructuredGrid::SafeDownCast(dataReader->GetOutput());
232 }
233 else
234 {
235 // for new filetypes
236 vtkXMLUnstructuredGridReader* xmlReader =
237 vtkXMLUnstructuredGridReader::SafeDownCast(algorithm);
238 grid = vtkUnstructuredGrid::SafeDownCast(xmlReader->GetOutput());
239 }
240 }
241 MeshLib::Mesh* mesh =
243 mesh->setName(item->data(0).toString().toStdString());
244 emit meshAdded(mesh);
245}
246
247void VtkVisPipelineView::selectionChanged(const QItemSelection& selected,
248 const QItemSelection& deselected)
249{
250 QTreeView::selectionChanged(selected, deselected);
251
252 if (selected.empty())
253 {
254 return;
255 }
256
257 QModelIndex index = *selected.indexes().begin();
258 if (index.isValid())
259 {
260 auto* item = static_cast<VtkVisPipelineItem*>(index.internalPointer());
261 emit actorSelected(item->actor());
262 emit itemSelected(item);
263 if (item->transformFilter())
264 {
265 emit dataObjectSelected(vtkDataObject::SafeDownCast(
266 item->transformFilter()->GetOutputDataObject(0)));
267 }
268 }
269 else
270 {
271 emit actorSelected(nullptr);
272 emit itemSelected(nullptr);
273 emit dataObjectSelected(nullptr);
274 }
275}
276
277void VtkVisPipelineView::selectItem(vtkProp3D* actor)
278{
279 this->selectItem(
280 static_cast<VtkVisPipeline*>(this->model())->getIndex(actor));
281}
282
283void VtkVisPipelineView::selectItem(const QModelIndex& index)
284{
285 if (!index.isValid())
286 {
287 return;
288 }
289
290 QItemSelectionModel* selectionModel = this->selectionModel();
291 selectionModel->clearSelection();
292 selectionModel->select(index, QItemSelectionModel::Select);
293}
294
296{
297 VtkVisPipelineItem* item(static_cast<VtkVisPipelineItem*>(
298 static_cast<VtkVisPipeline*>(this->model())
299 ->getItem(this->selectionModel()->currentIndex())));
300 const QString array_name = item->GetActiveAttribute();
301
302 QSettings settings;
303 QString filename = QFileDialog::getOpenFileName(
304 this, "Select color table",
305 settings.value("lastOpenedLutFileDirectory").toString(),
306 "Color table files (*.xml);;");
307 QFileInfo fi(filename);
308
309 if (fi.suffix().toLower() == "xml")
310 {
311 auto* pointSetItem = dynamic_cast<VtkVisPointSetItem*>(item);
312 if (pointSetItem)
313 {
314 VtkAlgorithmProperties* props = pointSetItem->getVtkProperties();
315 if (props)
316 {
317 props->SetLookUpTable(array_name, filename);
318 item->SetActiveAttribute(array_name);
319 emit requestViewUpdate();
320 }
321 }
322 else
323 {
324 QMessageBox::warning(nullptr,
325 "Color lookup table could not be applied.",
326 "Color lookup tables can only be applied to "
327 "VtkVisPointSetItem.");
328 }
329 QDir dir = QDir(filename);
330 settings.setValue("lastOpenedLutFileDirectory", dir.absolutePath());
331 }
332}
Definition of the CheckboxDelegate class.
Definition of the MeshFromRasterDialog class.
Definition of the Mesh class.
Definition of the OGSError class.
Definition of the VtkGeoImageSource class.
Definition of the VtkMeshConverter class.
Definition of the VtkVisImageItem class.
Definition of the VtkVisPipelineItem class.
Definition of the VtkVisPipelineView class.
Definition of the VtkVisPipeline class.
Definition of the VtkVisPointSetItem class.
CheckboxDelegate modifies a model view to display boolean values as checkboxes.
A dialog for specifying the parameters to construct a mesh based on a raster.
std::string getMeshName() const
MeshLib::MeshElemType getElementSelection() const
std::string getArrayName() const
MeshLib::UseIntensityAs getIntensitySelection() const
void setName(const std::string &name)
Changes the name of the mesh.
Definition Mesh.h:118
static MeshLib::Mesh * convertUnstructuredGrid(vtkUnstructuredGrid *grid, bool const compute_element_neighbors=false, std::string const &mesh_name="vtkUnstructuredGrid")
Converts a vtkUnstructuredGrid object to a Mesh.
static std::unique_ptr< MeshLib::Mesh > convert(GeoLib::Raster const &raster, MeshLib::MeshElemType elem_type, MeshLib::UseIntensityAs intensity_type, std::string const &array_name="Colour")
static void box(const QString &e)
Definition OGSError.cpp:23
Contains properties for the visualization of objects as VtkVisPipelineItems.
void SetLookUpTable(const QString &array_name, vtkLookupTable *lut)
Sets a colour lookup table for the given scalar array of the VtkVisPipelineItem.
bool IsRemovable() const
Is this algorithm removable from the pipeline (view).
static std::optional< GeoLib::Raster > convertToRaster(VtkGeoImageSource *const source)
An item in the VtkVisPipeline containing an image to be visualized.
An item in the VtkVisPipeline containing a graphic object to be visualized.
vtkAlgorithm * algorithm() const
Returns the algorithm object.
VtkAlgorithmProperties * getVtkProperties() const
Returns the VtkAlgorithmProperties.
virtual QString GetActiveAttribute() const
QVariant data(int column) const override
virtual void SetActiveAttribute(const QString &str)
void setModel(QAbstractItemModel *model) override
Overridden to set model specific header properties.
VtkVisPipelineView(QWidget *parent=nullptr)
Constructor.
void requestRemovePipelineItem(QModelIndex)
void meshAdded(MeshLib::Mesh *)
void dataObjectSelected(vtkDataObject *)
void requestAddPipelineFilterItem(QModelIndex)
void exportSelectedPipelineItemAsVtk()
Exports the currently selected item as a VTK file.
void convertVTKToOGSMesh()
Calls the conversion method for making a vtk grid an ogs mesh.
void actorSelected(vtkProp3D *)
void selectItem(vtkProp3D *actor)
void addColorTable()
Adds a color lookup table to the current scalar array of the selected pipeline item.
void itemSelected(VtkVisPipelineItem *)
void contextMenuEvent(QContextMenuEvent *event) override
Creates a menu on right-clicking on an item.
void showImageToMeshConversionDialog()
Calls the dialog to.
void writeRaster()
Calls the conversion method for saving this as an *.asc-file.
void addPipelineFilterItem()
Sends a requestAddPipelineFilterItem() signal to add a filter.
void selectionChanged(const QItemSelection &selected, const QItemSelection &deselected) override
Emits itemSelected() signals when an items was selected.
VtkVisPipeline manages the VTK visualization. It is a TreeModel and provides functions for adding and...
An item in the VtkVisPipeline containing a point set object to be visualized.
void writeToFile(std::string const &id_area_fname, std::string const &csv_fname, std::vector< std::pair< std::size_t, double > > const &ids_and_areas, std::vector< MeshLib::Node * > const &mesh_nodes)