ITK  6.0.0
Insight Toolkit
SphinxExamples/src/Core/ImageAdaptors/ProcessNthComponentOfVectorImage/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.
*
*=========================================================================*/
using VectorImageType = itk::Image<itk::CovariantVector<float, 3>, 2>;
static void
CreateImage(VectorImageType::Pointer image);
int
main()
{
auto image = VectorImageType::New();
CreateImage(image);
auto adaptor = ImageAdaptorType::New();
adaptor->SelectNthElement(0);
adaptor->SetImage(image);
auto blurFilter = BlurFilterType::New();
blurFilter->SetInput(adaptor);
blurFilter->Update();
return EXIT_SUCCESS;
}
void
CreateImage(VectorImageType::Pointer image)
{
auto size = VectorImageType::SizeType::Filled(2);
region.SetSize(size);
region.SetIndex(start);
image->SetRegions(region);
image->Allocate();
itk::ImageRegionIterator<VectorImageType> imageIterator(image, image->GetLargestPossibleRegion());
vec[0] = 1;
vec[1] = 2;
vec[2] = 3;
while (!imageIterator.IsAtEnd())
{
imageIterator.Set(vec);
++imageIterator;
}
}
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::NthElementImageAdaptor
Presents an image as being composed of the N-th element of its pixels.
Definition: itkNthElementImageAdaptor.h:58
itk::BinomialBlurImageFilter
Performs a separable blur on each dimension of an image.
Definition: itkBinomialBlurImageFilter.h:44
itkNthElementImageAdaptor.h
itkImageRegionIterator.h
itkImageAdaptor.h
itkBinomialBlurImageFilter.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::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
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