Skip to content
Snippets Groups Projects
Commit 455ce318 authored by Dhanannjay Deo's avatar Dhanannjay Deo Committed by Dhanannjay Deo
Browse files

FIX: Image streaming

- Release openslide handle correctly
- Also tests partial updates
- OpenSlideReader extent orientation
- Also uses new UpdateExtent to set reqeuested extents
parent c0615a90
No related branches found
No related tags found
No related merge requests found
# OpenSlideReader Tests
vtk_add_test_cxx(${vtk-module}CxxTests tests
TestOpenSlideReader.cxx
"DATA{${VTK_TEST_INPUT_DIR}/Microscopy/small2.ndpi}"
)
vtk_add_test_cxx(${vtk-module}CxxTests tests
TestOpenSlideReaderPartial.cxx
"DATA{${VTK_TEST_INPUT_DIR}/Microscopy/small2.ndpi}"
)
vtk_test_cxx_executable(${vtk-module}CxxTests tests RENDERING_FACTORY)
......@@ -41,15 +41,14 @@ int TestOpenSlideReader(int argc, char** argv)
reader->UpdateInformation();
delete [] rasterFileName;
int extent[6] = {20,120,20,120,0,0};
reader->SetUpdateExtent(extent);
vtkNew<vtkPNGWriter> writer;
writer->SetInputConnection(reader->GetOutputPort());
writer->SetFileName("this.png");
writer->SetUpdateExtent(extent);
writer->Update();
writer->Write();
// For debug
// reader->SetUpdateExtent(extent);
// vtkNew<vtkPNGWriter> writer;
// writer->SetInputConnection(reader->GetOutputPort());
// writer->SetFileName("this.png");
// writer->SetUpdateExtent(extent);
// writer->Update();
// writer->Write();
// Visualize
vtkNew<vtkRenderer> renderer;
......
/*=========================================================================
Program: Visualization Toolkit
Module: TestOpenSlideReader.cxx
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
#include <vtkNew.h>
#include <vtkOpenSlideReader.h>
#include <vtkRenderer.h>
#include <vtkRenderWindow.h>
#include <vtkRenderWindowInteractor.h>
#include <vtkImageViewer2.h>
#include <vtkImageData.h>
#include <vtkPNGWriter.h>
// VTK includes
#include <vtkTestUtilities.h>
// C++ includes
#include <sstream>
// Main program
int TestOpenSlideReaderPartial(int argc, char** argv)
{
const char* rasterFileName = vtkTestUtilities::ExpandDataFileName(argc, argv,
"Data/Microscopy/small2.ndpi");
//std::cout << "Got Filename: " << rasterFileName << std::endl;
// Create reader to read shape file.
vtkNew<vtkOpenSlideReader> reader;
reader->SetFileName(rasterFileName);
reader->UpdateInformation();
// reader->Print(cout);
delete [] rasterFileName;
int extent[6] = {200,499,200,499,0,0};
reader->UpdateExtent(extent);
vtkNew<vtkImageData> data;
data->ShallowCopy(reader->GetOutput());
// // For debug
// vtkNew<vtkPNGWriter> writer;
// writer->SetInputData(data.GetPointer());
// writer->SetFileName("this.png");
// writer->SetUpdateExtent(extent);
// writer->Update();
// writer->Write();
// Visualize
vtkNew<vtkRenderer> renderer;
vtkNew<vtkRenderWindow> window;
window->AddRenderer(renderer.GetPointer());
vtkNew<vtkRenderWindowInteractor> renderWindowInteractor;
renderWindowInteractor->SetRenderWindow(window.GetPointer());
vtkNew<vtkImageViewer2> imageViewer;
imageViewer->SetInputData(data.GetPointer());
//imageViewer->SetExtent(1000,1500,1000,1500,0,0);
imageViewer->SetupInteractor(renderWindowInteractor.GetPointer());
//imageViewer->SetSlice(0);
imageViewer->Render();
imageViewer->GetRenderer()->ResetCamera();
renderWindowInteractor->Initialize();
imageViewer->Render();
renderWindowInteractor->Start();
return EXIT_SUCCESS;
}
4b2972c7011281baa9a0fb3c883ffb7a
......@@ -85,10 +85,8 @@ void vtkOpenSlideReader::ExecuteDataWithInformation(vtkDataObject *output,
data->GetPointData()->GetScalars()->SetName("OpenSlideImage");
//// Leverage openslide to read the region
// int inExtent[6];
// data->GetExtent(inExtent);
// cout << inExtent[0] << ", " << inExtent[1] << endl;
// cout << inExtent[2] << ", " << inExtent[3] << endl;
// VTK extents have origin at top left and y axis looking downwards
// openslide needs to convert
int w = inExtent[1] - inExtent[0] + 1;
int h = inExtent[3]- inExtent[2] + 1;
......@@ -96,7 +94,7 @@ void vtkOpenSlideReader::ExecuteDataWithInformation(vtkDataObject *output,
openslide_read_region(this->openslide_handle, (unsigned int *) buffer,
inExtent[0],
inExtent[2],
this->DataExtent[3]-inExtent[3],
0, // level
w,
h
......@@ -104,7 +102,8 @@ void vtkOpenSlideReader::ExecuteDataWithInformation(vtkDataObject *output,
if(openslide_get_error(this->openslide_handle) != NULL)
{
delete[] buffer;
// Buffer is deleted by the openslide in case the error occurs
// delete[] buffer;
vtkErrorWithObjectMacro(this,
"File could not be read by openslide"
);
......@@ -129,7 +128,7 @@ void vtkOpenSlideReader::ExecuteDataWithInformation(vtkDataObject *output,
}
delete[] buffer;
openslide_close(this->openslide_handle);
// openslide_close(this->openslide_handle);
}
......@@ -150,10 +149,23 @@ int vtkOpenSlideReader::CanReadFile(const char* fname)
else
{
// Pretty sure
if(this->openslide_handle != NULL)
{
openslide_close(this->openslide_handle);
this->openslide_handle = NULL;
}
return 2;
}
}
vtkOpenSlideReader::~vtkOpenSlideReader()
{
// Release openslide_handle if being used
if(this->openslide_handle != NULL)
{
openslide_close(this->openslide_handle);
}
}
//----------------------------------------------------------------------------
void vtkOpenSlideReader::PrintSelf(ostream& os, vtkIndent indent)
......
......@@ -60,7 +60,7 @@ public:
}
protected:
vtkOpenSlideReader() {}
~vtkOpenSlideReader() {}
~vtkOpenSlideReader();
virtual void ExecuteInformation();
virtual void ExecuteDataWithInformation(vtkDataObject *out, vtkInformation *outInfo);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment