Add Constant to Every Pixel Without Duplicating Memory#
Synopsis#
Add a constant to every pixel without duplicating the image in memory (an accessor).
Results#
Output:
addPixelAccessor.SetValue(5)
image->GetPixel[0, 0]: 20 adaptor->GetPixel[0, 0]: 25
addPixelAccessor.SetValue(100)
image->GetPixel[0, 0]: 20 adaptor->GetPixel[0, 0]: 120
Code#
C++#
#include "itkImage.h"
#include "itkAddPixelAccessor.h"
#include "itkImageAdaptor.h"
#include "itkImageRegionIterator.h"
using ImageType = itk::Image<unsigned int, 2>;
static void
CreateImage(ImageType::Pointer image);
int
main()
{
auto image = ImageType::New();
CreateImage(image);
using AddPixelAccessorType = itk::Accessor::AddPixelAccessor<ImageType::PixelType>;
using ImageAdaptorType = itk::ImageAdaptor<ImageType, AddPixelAccessorType>;
auto adaptor = ImageAdaptorType::New();
AddPixelAccessorType addPixelAccessor;
adaptor->SetImage(image);
ImageType::IndexType index;
index[0] = 0;
index[1] = 0;
addPixelAccessor.SetValue(5);
adaptor->SetPixelAccessor(addPixelAccessor);
std::cout << "addPixelAccessor.SetValue(5)" << std::endl;
std::cout << "\timage->GetPixel" << index << ": " << image->GetPixel(index) << " adaptor->GetPixel" << index << ": "
<< adaptor->GetPixel(index) << std::endl;
addPixelAccessor.SetValue(100);
adaptor->SetPixelAccessor(addPixelAccessor);
std::cout << "addPixelAccessor.SetValue(100)" << std::endl;
std::cout << "\timage->GetPixel" << index << ": " << image->GetPixel(index) << " adaptor->GetPixel" << index << ": "
<< adaptor->GetPixel(index) << std::endl;
return EXIT_SUCCESS;
}
void
CreateImage(ImageType::Pointer image)
{
ImageType::IndexType start;
start.Fill(0);
ImageType::SizeType size;
size.Fill(10);
ImageType::RegionType region;
region.SetSize(size);
region.SetIndex(start);
image->SetRegions(region);
image->Allocate();
itk::ImageRegionIterator<ImageType> imageIterator(image, image->GetLargestPossibleRegion());
while (!imageIterator.IsAtEnd())
{
imageIterator.Set(20);
++imageIterator;
}
}
Classes demonstrated#
-
template<typename TPixel>
class AddPixelAccessor Simulates the effect of adding a constant value to all pixels.
This class is intended to be used as parameter of an ImageAdaptor to make an image appear as having pixels with intensity values increased by a constant amount.
- See
ImageAdaptor
- ITK Sphinx Examples: