ITK  6.0.0
Insight Toolkit
SphinxExamples/src/Core/ImageAdaptors/ViewComponentVectorImageAsScaleImage/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 "itkVectorImage.h"
#include "itkImage.h"
using ScalarImageType = itk::Image<float, 2>;
using VectorImageType = itk::VectorImage<float, 2>;
void
CreateImage(VectorImageType::Pointer image);
int
main()
{
auto image = VectorImageType::New();
CreateImage(image);
using ImageAdaptorType = itk::VectorImageToImageAdaptor<float, 2>;
auto adaptor = ImageAdaptorType::New();
adaptor->SetExtractComponentIndex(0);
adaptor->SetImage(image);
index[0] = 0;
index[1] = 0;
std::cout << adaptor->GetPixel(index) << std::endl;
return EXIT_SUCCESS;
}
void
CreateImage(VectorImageType::Pointer image)
{
start.Fill(0);
size.Fill(2);
region.SetSize(size);
region.SetIndex(start);
image->SetRegions(region);
image->SetNumberOfComponentsPerPixel(2);
image->Allocate();
using VariableVectorType = itk::VariableLengthVector<double>;
VariableVectorType variableLengthVector;
variableLengthVector.SetSize(2);
itk::ImageRegionIterator<VectorImageType> imageIterator(image, image->GetLargestPossibleRegion());
variableLengthVector[0] = 1;
variableLengthVector[1] = 2;
while (!imageIterator.IsAtEnd())
{
imageIterator.Set(variableLengthVector);
++imageIterator;
}
}
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itkVectorImageToImageAdaptor.h
itk::Index
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:68
itk::VectorImageToImageAdaptor
Presents a VectorImage and extracts a component from it into an image.
Definition: itkVectorImageToImageAdaptor.h:148
itk::VectorImage
Templated n-dimensional vector image class.
Definition: itkImageAlgorithm.h:29
itk::VariableLengthVector::SetSize
void SetSize(unsigned int sz, TReallocatePolicy reallocatePolicy, TKeepValuesPolicy keepValues)
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itkImageRegionIterator.h
itk::Index::Fill
void Fill(IndexValueType value)
Definition: itkIndex.h:272
itk::Size::Fill
void Fill(SizeValueType value)
Definition: itkSize.h:211
itkVectorImage.h
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::VariableLengthVector
Represents an array whose length can be defined at run-time.
Definition: itkConstantBoundaryCondition.h:28
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itk::ImageRegion::SetSize
void SetSize(const SizeType &size)
Definition: itkImageRegion.h:202