Convert Real and Imaginary Images to Complex Image#
Synopsis#
Convert a real image and an imaginary image to a complex image.
Results#
Output:
Image (0x7ff58250cf20)
RTTI typeinfo: itk::Image<std::__1::complex<float>, 2u>
Reference Count: 1
Modified Time: 25
Debug: Off
Object Name:
Observers:
none
Source: (0x7ff58250c590)
Source output name: Primary
Release Data: Off
Data Released: False
Global Release Data: Off
PipelineMTime: 22
UpdateMTime: 26
RealTimeStamp: 0 seconds
LargestPossibleRegion:
Dimension: 2
Index: [0, 0]
Size: [0, 0]
BufferedRegion:
Dimension: 2
Index: [0, 0]
Size: [0, 0]
RequestedRegion:
Dimension: 2
Index: [0, 0]
Size: [0, 0]
Spacing: [1, 1]
Origin: [0, 0]
Direction:
1 0
0 1
IndexToPointMatrix:
1 0
0 1
PointToIndexMatrix:
1 0
0 1
Inverse Direction:
1 0
0 1
PixelContainer:
ImportImageContainer (0x7ff58250d0e0)
RTTI typeinfo: itk::ImportImageContainer<unsigned long, std::__1::complex<float> >
Reference Count: 1
Modified Time: 24
Debug: Off
Object Name:
Observers:
none
Pointer: 0x7ff58250c1b0
Container manages memory: true
Size: 0
Capacity: 0
Code#
C++#
#include "itkImage.h"
#include "itkComposeImageFilter.h"
#include <complex>
int
main(int itkNotUsed(argc), char * itkNotUsed(argv)[])
{
using ImageType = itk::Image<unsigned char, 2>;
using ComplexImageType = itk::Image<std::complex<float>, 2>;
auto realImage = ImageType::New();
auto imaginaryImage = ImageType::New();
using RealAndImaginaryToComplexImageFilterType = itk::ComposeImageFilter<ImageType, ComplexImageType>;
RealAndImaginaryToComplexImageFilterType::Pointer realAndImaginaryToComplexImageFilter =
RealAndImaginaryToComplexImageFilterType::New();
realAndImaginaryToComplexImageFilter->SetInput1(realImage);
realAndImaginaryToComplexImageFilter->SetInput2(imaginaryImage);
realAndImaginaryToComplexImageFilter->Update();
ComplexImageType * output = realAndImaginaryToComplexImageFilter->GetOutput();
output->Print(std::cout);
return EXIT_SUCCESS;
}
Classes demonstrated#
-
template<typename TInputImage, typename TOutputImage = VectorImage<typename TInputImage::PixelType, TInputImage::ImageDimension>>
class ComposeImageFilter : public itk::ImageToImageFilter<TInputImage, TOutputImage> ComposeImageFilter combine several scalar images into a multicomponent image.
ComposeImageFilter combine several scalar images into an itk::Image of vector pixel (itk::Vector, itk::RGBPixel, …), of std::complex pixel, or in an itk::VectorImage.
- Inputs and Usage
- All input images are expected to have the same template parameters and have the same size and origin.
filter->SetInput( 0, image0 ); filter->SetInput( 1, image1 ); ... filter->Update(); itk::VectorImage< PixelType, dimension >::Pointer = filter->GetOutput();
- See
VectorImage
- See
VectorIndexSelectionCastImageFilter
- ITK Sphinx Examples: