Mask Image¶
Synopsis¶
Apply a mask to an image.
Results¶
Input Image¶
Output In VTK Window¶
Output:
Image (0x7f9d9c95ccc0)
RTTI typeinfo:   itk::Image<unsigned char, 2u>
Reference Count: 2
Modified Time: 235
Debug: Off
Object Name:
Observers:
  none
Source: (none)
Source output name: (none)
Release Data: Off
Data Released: False
Global Release Data: Off
PipelineMTime: 0
UpdateMTime: 0
RealTimeStamp: 0 seconds
LargestPossibleRegion:
  Dimension: 2
  Index: [0, 0]
  Size: [512, 342]
BufferedRegion:
  Dimension: 2
  Index: [0, 0]
  Size: [512, 342]
RequestedRegion:
  Dimension: 2
  Index: [0, 0]
  Size: [512, 342]
Spacing: [1, 1]
Origin: [0, 0]
Direction:
1 0
0 1
IndexToPointMatrix:
1 0
0 1
PointToIndexMatrix:
1 0
0 1
Inverse Direction:
1 0
0 1
PixelContainer:
  ImportImageContainer (0x7f9d9c95d000)
    RTTI typeinfo:   itk::ImportImageContainer<unsigned long, unsigned char>
    Reference Count: 1
    Modified Time: 236
    Debug: Off
    Object Name:
    Observers:
      none
    Pointer: 0x7f9d90050000
    Container manages memory: true
    Size: 175104
    Capacity: 175104
Code¶
C++¶
#include "itkConfigure.h"
#if (ITK_VERSION_MAJOR < 4) // These are all defaults in ITKv4
//  Not supported in ITKv3.
int
main(int argc, char * argv[])
{
  return 0;
}
#else
#  include "itkImage.h"
#  include "itkImageFileReader.h"
#  include "itkMaskImageFilter.h"
#  include "itkImageRegionIterator.h"
#  ifdef ENABLE_QUICKVIEW
#    include "QuickView.h"
#  endif
using ImageType = itk::Image<unsigned char, 2>;
void
CreateHalfMask(ImageType::Pointer image, ImageType::Pointer & mask);
int
main(int argc, char * argv[])
{
  if (argc < 2)
  {
    std::cerr << "Usage: " << argv[0] << " filename" << std::endl;
    return EXIT_FAILURE;
  }
  using ReaderType = itk::ImageFileReader<ImageType>;
  ReaderType::Pointer reader = ReaderType::New();
  reader->SetFileName(argv[1]);
  reader->Update();
  ImageType::Pointer mask = ImageType::New();
  CreateHalfMask(reader->GetOutput(), mask);
  using MaskFilterType = itk::MaskImageFilter<ImageType, ImageType>;
  MaskFilterType::Pointer maskFilter = MaskFilterType::New();
  maskFilter->SetInput(reader->GetOutput());
  maskFilter->SetMaskImage(mask);
  mask->Print(std::cout);
#  ifdef ENABLE_QUICKVIEW
  QuickView viewer;
  viewer.AddImage(reader->GetOutput(), true, itksys::SystemTools::GetFilenameName(argv[1]));
  std::stringstream desc;
  desc << "Mask";
  viewer.AddImage(mask.GetPointer(), true, desc.str());
  std::stringstream desc2;
  desc2 << "MaskFilter";
  viewer.AddImage(maskFilter->GetOutput(), true, desc2.str());
  viewer.Visualize();
#  endif
  return EXIT_SUCCESS;
}
void
CreateHalfMask(ImageType::Pointer image, ImageType::Pointer & mask)
{
  ImageType::RegionType region = image->GetLargestPossibleRegion();
  mask->SetRegions(region);
  mask->Allocate();
  ImageType::SizeType regionSize = region.GetSize();
  itk::ImageRegionIterator<ImageType> imageIterator(mask, region);
  // Make the left half of the mask white and the right half black
  while (!imageIterator.IsAtEnd())
  {
    if (imageIterator.GetIndex()[0] > static_cast<ImageType::IndexValueType>(regionSize[0]) / 2)
    {
      imageIterator.Set(0);
    }
    else
    {
      imageIterator.Set(255);
    }
    ++imageIterator;
  }
}
#endif
Classes demonstrated¶
- 
template<typename 
TInputImage, typenameTMaskImage, typenameTOutputImage= TInputImage>
classMaskImageFilter: public itk::BinaryGeneratorImageFilter<TInputImage, TMaskImage, TOutputImage> Mask an image with a mask.
This class is templated over the types of the input image type, the mask image type and the type of the output image. Numeric conversions (castings) are done by the C++ defaults.
The pixel type of the input 2 image must have a valid definition of the operator != with zero. This condition is required because internally this filter will perform the operation
if pixel_from_mask_image != masking_value pixel_output_image = pixel_input_image else pixel_output_image = outside_value
The pixel from the input 1 is cast to the pixel type of the output image.
Note that the input and the mask images must be of the same size.
- Warning
 Any pixel value other than masking value (0 by default) will not be masked out.
- See
 MaskNegatedImageFilter
- ITK Sphinx Examples:
 

