Stream a Pipeline#
Synopsis#
Stream a pipeline. The PipelineMonitorImageFilter shows how the image is processed in three sub-regions. Note that the StreamingImageFilter comes at the end of the pipeline and that no calls to Update() occur on any intermediate filters.
Results#
The output LargestPossibleRegion is: ImageRegion (0x8dc420)
Dimension: 2
Index: [0, 0]
Size: [3, 3]
Updated RequestedRegions's:
ImageRegion (0x8dce80)
Dimension: 2
Index: [0, 0]
Size: [3, 1]
ImageRegion (0x8dcea8)
Dimension: 2
Index: [0, 1]
Size: [3, 1]
ImageRegion (0x8dced0)
Dimension: 2
Index: [0, 2]
Size: [3, 1]
Code#
Python#
#!/usr/bin/env python
import sys
import itk
import argparse
from distutils.version import StrictVersion as VS
if VS(itk.Version.GetITKVersion()) < VS("4.10.0"):
print("ITK 4.10.0 is required.")
sys.exit(1)
parser = argparse.ArgumentParser(description="Stream A Pipeline.")
parser.add_argument("number_of_splits", type=int)
args = parser.parse_args()
Dimension = 2
PixelType = itk.UC
ImageType = itk.Image[PixelType, Dimension]
source = itk.RandomImageSource[ImageType].New()
size = itk.Size[Dimension]()
size.Fill(args.number_of_splits)
source.SetSize(size)
monitorFilter = itk.PipelineMonitorImageFilter[ImageType].New()
monitorFilter.SetInput(source.GetOutput())
streamingFilter = itk.StreamingImageFilter[ImageType, ImageType].New()
streamingFilter.SetInput(monitorFilter.GetOutput())
streamingFilter.SetNumberOfStreamDivisions(args.number_of_splits)
streamingFilter.Update()
print(
"The output LargestPossibleRegion is: "
+ str(streamingFilter.GetOutput().GetLargestPossibleRegion())
)
print("")
updatedRequestedRegions = monitorFilter.GetUpdatedRequestedRegions()
print("Updated ApplyAFilterOnlyToASpecifiedImageRegion's:")
for ii in range(len(updatedRequestedRegions)):
print(" " + str(updatedRequestedRegions[ii]))
C++#
#include "itkPipelineMonitorImageFilter.h"
#include "itkRandomImageSource.h"
#include "itkStreamingImageFilter.h"
int
main(int argc, char * argv[])
{
if (argc != 2)
{
std::cerr << "Usage: " << argv[0] << " <NumberOfSplits>" << std::endl;
return EXIT_FAILURE;
}
int numberOfSplits = std::stoi(argv[1]);
constexpr unsigned int Dimension = 2;
using PixelType = unsigned char;
using ImageType = itk::Image<PixelType, Dimension>;
using SourceType = itk::RandomImageSource<ImageType>;
auto source = SourceType::New();
ImageType::SizeType size;
size.Fill(numberOfSplits);
source->SetSize(size);
using MonitorFilterType = itk::PipelineMonitorImageFilter<ImageType>;
auto monitorFilter = MonitorFilterType::New();
monitorFilter->SetInput(source->GetOutput());
// If ITK was built with the Debug CMake configuration, the filter
// automatically outputs status information to the console
monitorFilter->DebugOn();
using StreamingFilterType = itk::StreamingImageFilter<ImageType, ImageType>;
auto streamingFilter = StreamingFilterType::New();
streamingFilter->SetInput(monitorFilter->GetOutput());
streamingFilter->SetNumberOfStreamDivisions(numberOfSplits);
try
{
streamingFilter->Update();
}
catch (const itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
std::cout << "The output LargestPossibleRegion is: " << streamingFilter->GetOutput()->GetLargestPossibleRegion()
<< std::endl;
std::cout << std::endl;
const MonitorFilterType::RegionVectorType updatedRequestedRegions = monitorFilter->GetUpdatedRequestedRegions();
std::cout << "Updated RequestedRegions's:" << std::endl;
for (const auto & updatedRequestedRegion : updatedRequestedRegions)
{
std::cout << " " << updatedRequestedRegion << std::endl;
}
return EXIT_SUCCESS;
}
Classes demonstrated#
-
template<typename TInputImage, typename TOutputImage>
class StreamingImageFilter : public itk::ImageToImageFilter<TInputImage, TOutputImage> Pipeline object to control data streaming for large data processing.
StreamingImageFilter is a pipeline object that allows the user to control how data is pulled through the pipeline. To generate its OutputRequestedRegion, this filter will divide the output into several pieces (controlled by SetNumberOfStreamDivisions), and call the upstream pipeline for each piece, tiling the individual outputs into one large output. This reduces the memory footprint for the application since each filter does not have to process the entire dataset at once. This filter will produce the entire output as one image, but the upstream filters will do their processing in pieces.
-
template<typename TImageType>
class PipelineMonitorImageFilter : public itk::ImageToImageFilter<TImageType, TImageType> Enables monitoring, recording and debugging of the pipeline execution and information exchange.
This filter is useful for testing, debugging, and understanding the pipeline. When DebugOn is enabled and compiled in Debug mode, many itkDebug messages are printed. This filter also features, several Verify methods which check the recorded information, for certain conditions, which should occur when well behaved filters are executed.
There are two meta verify methods that should primarily be used depending on the expected capabilities of the pipeline:
VerifyAllInputCanStream
VerifyAllInputCanNotStream
During the pipeline execution this filter records a variety of information to aid if verifying correct pipeline behavior:
NumberOfUpdate the number of times GenerateData was executed
OutputRequestedRegions an array of the output of this filter’s requested region
InputRequestedRegions an array of the input image’s requested region after PropagateRequestedRegion
UpdatedBufferedRegions an array of the input image’s buffered region after upstream GenerateData
UpdatedRequestedRegions an array of the input image’s requested region after upstream GenerateData
The following are recorded from the input image after the input’s output information is generated:
UpdatedOutputOrigin
UpdatedOutputDirection
UpdatedOutputSpacing
UpdatedOutputLargestPossibleRegion
This filter always runs in-place so it has no per-pixel overhead.