Watch a Filter#

Synopsis#

This example demonstrates how to watch what is going on inside of a filter.

Results#

Code#

C++#

#include "itkImage.h"
#include "itkSimpleFilterWatcher.h"
#include "itkThresholdImageFilter.h"
#include "itkImageRegionIterator.h"

template <class TImage>
void
CreateImage(typename TImage::Pointer image)
{
  using ImageType = TImage;

  // Create an image with 2 connected components
  typename ImageType::RegionType region;
  typename ImageType::IndexType  start;
  start.Fill(0);

  typename ImageType::SizeType size;
  size.Fill(100);

  region.SetSize(size);
  region.SetIndex(start);

  image->SetRegions(region);
  image->Allocate();

  itk::ImageRegionIterator<ImageType> imageIterator(image, region);

  while (!imageIterator.IsAtEnd())
  {
    imageIterator.Set(255);
    ++imageIterator;
  }
}
int
main()
{
  constexpr unsigned int Dimension = 2;
  using PixelType = unsigned char;

  using ImageType = itk::Image<PixelType, Dimension>;

  auto image = ImageType::New();
  CreateImage<ImageType>(image);

  using ThresholdImageFilterType = itk::ThresholdImageFilter<ImageType>;

  auto thresholdFilter = ThresholdImageFilterType::New();
  thresholdFilter->SetInput(image);
  thresholdFilter->ThresholdBelow(100);
  thresholdFilter->SetOutsideValue(0);

  itk::SimpleFilterWatcher watcher(thresholdFilter, "ThresholdFilter");

  thresholdFilter->Update();

  return EXIT_SUCCESS;
}

Classes demonstrated#

class SimpleFilterWatcher

Simple mechanism for monitoring the pipeline events of a filter and reporting these events to std::cout.

SimpleFilterWatcher provides a simple mechanism for monitoring the execution of filter. SimpleFilterWatcher is a stack-based object which takes a pointer to a ProcessObject at constructor time. SimpleFilterWatcher creates a series of commands that are registered as observers to the specified ProcessObject. The events monitored are:

 StartEvent
 EndEvent
 ProgressEvent
 IterationEvent
 AbortEvent

The callbacks routines registered for these events emit a simple message to std::cout.

Example of use:

using FilterType = itk::BinaryThresholdImageFilter<ImageType>; FilterType::Pointer thresholdFilter = FilterType::New();

SimpleFilterWatcher watcher(thresholdFilter, “Threshold”);

The second argument to the constructor to SimpleFilterWatcher is an optional string that is prepended to the event messages. This allows the user to associate the emitted messages to a particular filter/variable.

ITK Sphinx Examples:

Subclassed by itk::XMLFilterWatcher

See itk::SimpleFilterWatcher for additional documentation.