ITK  5.4.0
Insight Toolkit
SphinxExamples/src/Core/ImageAdaptors/AddConstantToPixelsWithoutDuplicatingImage/Code.cxx
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* https://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#include "itkImage.h"
using ImageType = itk::Image<unsigned int, 2>;
static void
CreateImage(ImageType::Pointer image);
int
main()
{
auto image = ImageType::New();
CreateImage(image);
auto adaptor = ImageAdaptorType::New();
AddPixelAccessorType addPixelAccessor;
adaptor->SetImage(image);
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)
{
start.Fill(0);
size.Fill(10);
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;
}
}
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itkAddPixelAccessor.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itkImageRegionIterator.h
itkImageAdaptor.h
itk::Index::Fill
void Fill(IndexValueType value)
Definition: itkIndex.h:274
itk::Size::Fill
void Fill(SizeValueType value)
Definition: itkSize.h:213
itk::ImageAdaptor
Give access to partial aspects of voxels from an Image.
Definition: itkImageAdaptor.h:56
itk::ImageRegionIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionIterator.h:80
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ImageRegionIterator::Set
void Set(const PixelType &value) const
Definition: itkImageRegionIterator.h:116
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itk::Accessor::AddPixelAccessor
Simulates the effect of adding a constant value to all pixels.
Definition: itkAddPixelAccessor.h:45
itk::ImageRegion::SetSize
void SetSize(const SizeType &size)
Definition: itkImageRegion.h:202