OGS
VtkVisImageItem.cpp
Go to the documentation of this file.
1 
15 // ** INCLUDES **
16 #include "VtkVisImageItem.h"
17 
18 #include <vtkActor.h>
19 #include <vtkDataArray.h>
20 #include <vtkDataSetMapper.h>
21 #include <vtkImageAlgorithm.h>
22 #include <vtkImageChangeInformation.h>
23 #include <vtkImageData.h>
24 #include <vtkImageShiftScale.h>
25 #include <vtkPointData.h>
26 #include <vtkRenderer.h>
27 #include <vtkSmartPointer.h>
28 
29 #include "BaseLib/FileTools.h"
30 #include "BaseLib/Logging.h"
31 #include "VtkAlgorithmProperties.h"
32 #include "VtkGeoImageSource.h"
33 
34 // export
35 #include <vtkImageActor.h>
36 #include <vtkXMLImageDataWriter.h>
37 
39  vtkAlgorithm* algorithm, TreeItem* parentItem,
40  const QList<QVariant> data /*= QList<QVariant>()*/)
41  : VtkVisPipelineItem(algorithm, parentItem, data), _transformFilter(nullptr)
42 {
43 }
44 
46  VtkCompositeFilter* compositeFilter, TreeItem* parentItem,
47  const QList<QVariant> data /*= QList<QVariant>()*/)
48  : VtkVisPipelineItem(compositeFilter, parentItem, data),
49  _transformFilter(nullptr)
50 {
51 }
52 
54 {
55  _transformFilter->Delete();
56 }
57 
58 void VtkVisImageItem::Initialize(vtkRenderer* renderer)
59 {
60  auto* img = dynamic_cast<vtkImageAlgorithm*>(_algorithm);
61  img->Update();
62  // VtkGeoImageSource* img = dynamic_cast<VtkGeoImageSource*>(_algorithm);
63 
64  double origin[3];
65  double spacing[3];
66  double range[2];
67  img->GetOutput()->GetOrigin(origin);
68  img->GetOutput()->GetSpacing(spacing);
69  // img->getRange(range);
70  img->GetOutput()->GetPointData()->GetScalars()->GetRange(range);
71  vtkImageShiftScale* scale = vtkImageShiftScale::New();
72  scale->SetOutputScalarTypeToUnsignedChar();
73  scale->SetInputConnection(img->GetOutputPort());
74  scale->SetShift(-range[0]);
75  scale->SetScale(255.0 / (range[1] - range[0]));
76 
77  _transformFilter = vtkImageChangeInformation::New();
78  _transformFilter->SetInputConnection(scale->GetOutputPort());
79  // double origin[3];
80  // img->getOrigin(origin);
81  // double spacing = img->getSpacing();
82  //_transformFilter->SetOutputOrigin(origin2);
83  //_transformFilter->SetOutputSpacing(spacing2);
84  _transformFilter->Update();
85 
86  _renderer = renderer;
87 
88  // Use a special vtkImageActor instead of vtkActor
89  vtkImageActor* imageActor = vtkImageActor::New();
90  imageActor->SetInputData(_transformFilter->GetOutput());
91  _actor = imageActor;
92  _renderer->AddActor(_actor);
93 
94  // Set pre-set properties
95  auto* vtkProps = dynamic_cast<VtkAlgorithmProperties*>(_algorithm);
96  if (vtkProps)
97  {
98  setVtkProperties(vtkProps);
99  }
100 
101  auto* parentItem = dynamic_cast<VtkVisPipelineItem*>(this->parentItem());
102  while (parentItem)
103  {
104  auto* parentProps =
105  dynamic_cast<VtkAlgorithmProperties*>(parentItem->algorithm());
106  if (parentProps)
107  {
108  auto* newProps = new VtkAlgorithmProperties();
109  newProps->SetScalarVisibility(parentProps->GetScalarVisibility());
110  newProps->SetTexture(parentProps->GetTexture());
111  setVtkProperties(newProps);
112  vtkProps = newProps;
113  parentItem = nullptr;
114  }
115  else
116  {
117  parentItem =
118  dynamic_cast<VtkVisPipelineItem*>(parentItem->parentItem());
119  }
120  }
121 
122  // Set active scalar to the desired one from VtkAlgorithmProperties
123  // or to match those of the parent.
124  if (vtkProps)
125  {
126  if (vtkProps->GetActiveAttribute().length() > 0)
127  {
128  this->SetActiveAttribute(vtkProps->GetActiveAttribute());
129  }
130  else
131  {
132  auto* visParentItem =
133  dynamic_cast<VtkVisPipelineItem*>(this->parentItem());
134  if (visParentItem)
135  {
136  this->SetActiveAttribute(visParentItem->GetActiveAttribute());
137  }
138  if (vtkProps->GetTexture() != nullptr)
139  {
140  this->SetActiveAttribute("Solid Color");
141  }
142  }
143  }
144 }
145 
147 {
148  // todo implementation
149  (void)vtkProps;
150 }
151 
152 int VtkVisImageItem::callVTKWriter(vtkAlgorithm* algorithm,
153  const std::string& filename) const
154 {
155  std::string file_name_cpy(filename);
156  auto* algID = dynamic_cast<vtkImageAlgorithm*>(algorithm);
157  if (algID)
158  {
159  vtkSmartPointer<vtkXMLImageDataWriter> iWriter =
160  vtkSmartPointer<vtkXMLImageDataWriter>::New();
161  iWriter->SetInputData(algID->GetOutputDataObject(0));
162  if (BaseLib::getFileExtension(filename) != ".vti")
163  {
164  file_name_cpy.append(".vti");
165  }
166  iWriter->SetFileName(file_name_cpy.c_str());
167  return iWriter->Write();
168  }
169  ERR("VtkVisPipelineItem::writeToFile() - Unknown data type.");
170  return 0;
171 }
172 
173 void VtkVisImageItem::setTranslation(double x, double y, double z) const
174 {
175  _transformFilter->SetOriginTranslation(x, y, z);
176 }
177 
178 vtkAlgorithm* VtkVisImageItem::transformFilter() const
179 {
180  return _transformFilter;
181 }
Filename manipulation routines.
void ERR(char const *fmt, Args const &... args)
Definition: Logging.h:42
Definition of the VtkAlgorithmProperties class.
Definition of the VtkGeoImageSource class.
Definition of the VtkVisImageItem class.
Objects nodes for the TreeModel.
Definition: TreeItem.h:28
TreeItem * parentItem() const
Definition: TreeItem.cpp:115
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...
vtkImageChangeInformation * _transformFilter
VtkVisImageItem(vtkAlgorithm *algorithm, TreeItem *parentItem, const QList< QVariant > data=QList< QVariant >())
Constructor for a source/filter object.
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...
An item in the VtkVisPipeline containing a graphic object to be visualized.
vtkAlgorithm * algorithm() const
Returns the algorithm object.
vtkAlgorithm * _algorithm
virtual void SetActiveAttribute(const QString &str)
std::string getFileExtension(const std::string &path)
Definition: FileTools.cpp:186
void scale(PETScVector &x, double const a)
Definition: LinAlg.cpp:44