Observe an Event#

Synopsis#

This example demonstrates how to observe and event that is invoked by a filter.

Results#

Output:

Progress: 0
Progress: 0.00994873
Progress: 0.0198975
Progress: 0.0298462
Progress: 0.0397949
Progress: 0.0497437
Progress: 0.0596924
Progress: 0.0696411
Progress: 0.0795898
Progress: 0.0895386
Progress: 0.0994873
Progress: 0.109436
Progress: 0.119385
Progress: 0.129333
Progress: 0.139282
Progress: 0.149231
Progress: 0.15918
Progress: 0.169128
Progress: 0.179077
Progress: 0.189026
Progress: 0.198975
Progress: 0.208923
Progress: 0.218872
Progress: 0.228821
Progress: 0.23877
Progress: 0.248718
Progress: 0.258667
Progress: 0.268616
Progress: 0.278564
Progress: 0.288513
Progress: 0.298462
Progress: 0.308411
Progress: 0.318359
Progress: 0.328308
Progress: 0.338257
Progress: 0.348206
Progress: 0.358154
Progress: 0.368103
Progress: 0.378052
Progress: 0.388
Progress: 0.397949
Progress: 0.407898
Progress: 0.417847
Progress: 0.427795
Progress: 0.437744
Progress: 0.447693
Progress: 0.457642
Progress: 0.46759
Progress: 0.477539
Progress: 0.487488
Progress: 0.497437
Progress: 0.507385
Progress: 0.517334
Progress: 0.527283
Progress: 0.537231
Progress: 0.54718
Progress: 0.557129
Progress: 0.567078
Progress: 0.577026
Progress: 0.586975
Progress: 0.596924
Progress: 0.606873
Progress: 0.616821
Progress: 0.62677
Progress: 0.636719
Progress: 0.646667
Progress: 0.656616
Progress: 0.666565
Progress: 0.676514
Progress: 0.686462
Progress: 0.696411
Progress: 0.70636
Progress: 0.716309
Progress: 0.726257
Progress: 0.736206
Progress: 0.746155
Progress: 0.756104
Progress: 0.766052
Progress: 0.776001
Progress: 0.78595
Progress: 0.795898
Progress: 0.805847
Progress: 0.815796
Progress: 0.825745
Progress: 0.835693
Progress: 0.845642
Progress: 0.855591
Progress: 0.86554
Progress: 0.875488
Progress: 0.885437
Progress: 0.895386
Progress: 0.905334
Progress: 0.915283
Progress: 0.925232
Progress: 0.935181
Progress: 0.945129
Progress: 0.955078
Progress: 0.965027
Progress: 0.974976
Progress: 0.984924
Progress: 0.994873
Progress: 1

Code#

Python#

#!/usr/bin/env python

import itk

Dimension = 2
PixelType = itk.UC
ImageType = itk.Image[PixelType, Dimension]

source = itk.GaussianImageSource[ImageType].New()
size = itk.Size[Dimension]()
size.Fill(128)
source.SetSize(size)

sigma = itk.FixedArray[itk.D, Dimension]()
sigma.Fill(45.0)
source.SetSigma(sigma)


def myCommand():
    print("Progress: " + str(source.GetProgress()))


source.AddObserver(itk.ProgressEvent(), myCommand)

source.Update()

C++#

#include "itkImage.h"
#include "itkGaussianImageSource.h"
#include "itkCommand.h"

class MyCommand : public itk::Command
{
public:
  itkNewMacro(MyCommand);

  void
  Execute(itk::Object * caller, const itk::EventObject & event) override
  {
    Execute((const itk::Object *)caller, event);
  }

  void
  Execute(const itk::Object * caller, const itk::EventObject & event) override
  {
    if (!itk::ProgressEvent().CheckEvent(&event))
    {
      return;
    }
    const auto * processObject = dynamic_cast<const itk::ProcessObject *>(caller);
    if (!processObject)
    {
      return;
    }
    std::cout << "Progress: " << processObject->GetProgress() << std::endl;
  }
};


int
main()
{
  constexpr unsigned int Dimension = 2;
  using PixelType = unsigned char;
  using ImageType = itk::Image<PixelType, Dimension>;

  using SourceType = itk::GaussianImageSource<ImageType>;
  auto source = SourceType::New();

  ImageType::SizeType size;
  size.Fill(128);
  source->SetSize(size);

  SourceType::ArrayType sigma;
  sigma.Fill(45.0);
  source->SetSigma(sigma);

  auto myCommand = MyCommand::New();
  source->AddObserver(itk::ProgressEvent(), myCommand);

  source->Update();

  return EXIT_SUCCESS;
}

Classes demonstrated#

class Command : public itk::Object

Superclass for callback/observer methods.

Command is an implementation of the command design pattern that is used in callbacks (such as StartMethod(), ProgressMethod(), and EndMethod()) in ITK. itk::Object implements a Subject/Observer pattern. When a subject needs to notify a observer, it does so using a itk::Command. The Execute method is called to run the command.

ITK Sphinx Examples:

Subclassed by itk::CommandIterationUpdate< TOptimizer >, itk::CommandIterationUpdatev4< TOptimizer >, itk::CommandVnlIterationUpdate< TOptimizer >, itk::CStyleCommand, itk::FunctionCommand, itk::MemberCommand< T >, itk::ReceptorMemberCommand< T >, itk::SimpleConstMemberCommand< T >, itk::SimpleMemberCommand< T >, itk::WatershedMiniPipelineProgressCommand

See itk::Command for additional documentation.