ITK
6.0.0
Insight Toolkit
Examples/DataRepresentation/Image/ImageAdaptor3.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.
*
*=========================================================================*/
// Software Guide : BeginLatex
//
// This example illustrates the use of \doxygen{ImageAdaptor}
// to obtain access to the components of a vector image.
// Specifically, it shows how to manage pixel accessors containing
// internal parameters. In this example we create an image of vectors by using
// a gradient filter. Then, we use an image adaptor to extract one of the
// components of the vector image. The vector type used by the gradient filter
// is the \doxygen{CovariantVector} class.
//
// We start by including the relevant headers.
//
// \index{itk::ImageAdaptor!Instantiation}
// \index{itk::ImageAdaptor!Header}
// \index{itk::PixelAccessor!with parameters}
//
// Software Guide : EndLatex
#include "
itkImageAdaptor.h
"
#include "
itkImageRegionIteratorWithIndex.h
"
#include "
itkImageFileReader.h
"
#include "
itkImageFileWriter.h
"
#include "
itkRescaleIntensityImageFilter.h
"
// Software Guide : BeginCodeSnippet
#include "
itkGradientRecursiveGaussianImageFilter.h
"
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// A pixel accessors class may have internal parameters that affect the
// operations performed on input pixel data. Image adaptors support
// parameters in their internal pixel accessor by using
// the assignment operator. Any pixel accessor which has internal
// parameters must therefore implement the assignment operator.
// The following defines a pixel accessor for extracting
// components from a vector pixel. The
// \code{m\_Index} member variable is used to select the vector component
// to be returned.
//
// Software Guide : EndLatex
namespace
itk
{
// Software Guide : BeginCodeSnippet
class
VectorPixelAccessor
{
public
:
using
InternalType =
itk::CovariantVector<float, 2>
;
using
ExternalType = float;
VectorPixelAccessor() =
default
;
VectorPixelAccessor &
operator=(
const
VectorPixelAccessor & vpa) =
default
;
ExternalType
Get(
const
InternalType & input)
const
{
return
static_cast<ExternalType>(input[m_Index]);
}
void
SetIndex(
unsigned
int
index)
{
m_Index = index;
}
private
:
unsigned
int
m_Index{ 0 };
};
// Software Guide : EndCodeSnippet
}
// namespace itk
// Software Guide : BeginLatex
//
// The \code{Get()} method simply returns the \emph{i}-th component of
// the vector as indicated by the index. The assignment operator transfers
// the value of the index member variable from one instance of the pixel
// accessor to another.
//
// Software Guide : EndLatex
//-------------------------
//
// Main code
//
//-------------------------
int
main(
int
argc,
char
* argv[])
{
if
(argc < 4)
{
std::cerr <<
"Usage: "
<< std::endl;
std::cerr <<
"ImageAdaptor3 inputFileName outputComponentFileName "
;
std::cerr <<
" indexOfComponentToExtract"
<< std::endl;
return
EXIT_FAILURE;
}
// Software Guide : BeginLatex
//
// In order to test the pixel accessor, we generate an image of vectors
// using the \doxygen{GradientRecursiveGaussianImageFilter}. This filter
// produces an output image of \doxygen{CovariantVector} pixel type.
// Covariant vectors are the natural representation for gradients since
// they are the equivalent of normals to iso-values manifolds.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
using
InputPixelType =
unsigned
char;
constexpr
unsigned
int
Dimension
= 2;
using
InputImageType =
itk::Image<InputPixelType, Dimension>
;
using
VectorPixelType =
itk::CovariantVector<float, Dimension>
;
using
VectorImageType =
itk::Image<VectorPixelType, Dimension>
;
using
GradientFilterType =
itk::GradientRecursiveGaussianImageFilter
<InputImageType,
VectorImageType>;
auto
gradient =
GradientFilterType::New
();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We instantiate the ImageAdaptor using the vector image type as
// the first template parameter and the pixel accessor as the second
// template parameter.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
using
ImageAdaptorType =
itk::ImageAdaptor<VectorImageType, itk::VectorPixelAccessor>
;
auto
adaptor =
ImageAdaptorType::New
();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The index of the component to be extracted is specified
// from the command line. In the following, we create the accessor,
// set the index and connect the accessor to the image adaptor using
// the \code{SetPixelAccessor()} method.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
itk::VectorPixelAccessor accessor;
accessor.SetIndex(std::stoi(argv[3]));
adaptor->SetPixelAccessor(accessor);
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We create a reader to load the image specified from the
// command line and pass its output as the input to the gradient filter.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
using
ReaderType =
itk::ImageFileReader<InputImageType>
;
auto
reader =
ReaderType::New
();
gradient->SetInput(reader->GetOutput());
reader->SetFileName(argv[1]);
gradient->Update();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We now connect the output of the gradient filter as input to the
// image adaptor. The adaptor emulates a scalar image whose pixel values
// are taken from the selected component of the vector image.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
adaptor->SetImage(gradient->GetOutput());
// Software Guide : EndCodeSnippet
using
OutputImageType =
itk::Image<unsigned char, Dimension>
;
using
RescalerType =
itk::RescaleIntensityImageFilter<ImageAdaptorType, OutputImageType>
;
auto
rescaler =
RescalerType::New
();
using
WriterType =
itk::ImageFileWriter<OutputImageType>
;
auto
writer =
WriterType::New
();
writer->SetFileName(argv[2]);
rescaler->SetOutputMinimum(0);
rescaler->SetOutputMaximum(255);
rescaler->SetInput(adaptor);
writer->SetInput(rescaler->GetOutput());
writer->Update();
// Software Guide : BeginLatex
//
// \begin{figure} \center
// \includegraphics[width=0.32\textwidth]{BrainProtonDensitySlice}
// \includegraphics[width=0.32\textwidth]{ImageAdaptorToVectorImageComponentX}
// \includegraphics[width=0.32\textwidth]{ImageAdaptorToVectorImageComponentY}
// \itkcaption[Image Adaptor to Vector Image]{Using
// ImageAdaptor to access components of a vector
// image. The input image on the left was passed through a gradient image
// filter and the two components of the resulting vector image were
// extracted using an image adaptor.} \label{fig:ImageAdaptorToVectorImage}
// \end{figure}
//
// As in the previous example, we rescale the scalar image before writing
// the image out to file. Figure~\ref{fig:ImageAdaptorToVectorImage}
// shows the result of applying the example code for extracting both
// components of a two dimensional gradient.
//
// Software Guide : EndLatex
return
EXIT_SUCCESS;
}
itkImageFileReader.h
itkImageAdaptor.h
itkImageRegionIteratorWithIndex.h
itk::ImageAdaptor
Give access to partial aspects of voxels from an Image.
Definition:
itkImageAdaptor.h:56
itk::ImageFileReader
Data source that reads image data from a single file.
Definition:
itkImageFileReader.h:75
itkGradientRecursiveGaussianImageFilter.h
itk::ImageFileWriter
Writes image data to a single file.
Definition:
itkImageFileWriter.h:90
itk::GradientRecursiveGaussianImageFilter
Computes the gradient of an image by convolution with the first derivative of a Gaussian.
Definition:
itkGradientRecursiveGaussianImageFilter.h:59
itkRescaleIntensityImageFilter.h
itkImageFileWriter.h
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition:
itkCovariantVector.h:70
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition:
itkAnatomicalOrientation.h:29
itk::RescaleIntensityImageFilter
Applies a linear transformation to the intensity levels of the input Image.
Definition:
itkRescaleIntensityImageFilter.h:133
itk::Image
Templated n-dimensional image class.
Definition:
itkImage.h:88
New
static Pointer New()
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition:
itkGTestTypedefsAndConstructors.h:44
Generated on
unknown
for ITK by
1.8.16