ITK  6.0.0
Insight Toolkit
SphinxExamples/src/Core/Transform/ApplyAffineTransformFromHomogeneousMatrixAndResample/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.
*
*=========================================================================*/
int
main(int argc, char * argv[])
{
if (argc != 4)
{
std::cerr << "Usage: " << std::endl;
std::cerr << argv[0];
std::cerr << " <InputFileName> <OutputFileName> <DefaultPixelValue>";
std::cerr << std::endl;
return EXIT_FAILURE;
}
constexpr unsigned int Dimension = 2;
using ScalarType = double;
const char * inputFileName = argv[1];
const char * outputFileName = argv[2];
const auto defaultValue = static_cast<ScalarType>(std::stod(argv[3]));
MatrixType matrix;
matrix[0][0] = std::cos(0.05);
matrix[0][1] = std::sin(0.05);
matrix[0][2] = 0.;
matrix[1][0] = -matrix[0][1];
matrix[1][1] = matrix[0][0];
matrix[1][2] = 0.;
matrix[2][0] = -10.;
matrix[2][1] = -20.;
matrix[2][2] = 1.;
using PixelType = unsigned char;
const auto input = itk::ReadImage<ImageType>(inputFileName);
const ImageType::SizeType & size = input->GetLargestPossibleRegion().GetSize();
using ResampleImageFilterType = itk::ResampleImageFilter<ImageType, ImageType>;
auto resample = ResampleImageFilterType::New();
resample->SetInput(input);
resample->SetReferenceImage(input);
resample->UseReferenceImageOn();
resample->SetSize(size);
resample->SetDefaultPixelValue(defaultValue);
constexpr unsigned int Radius = 3;
auto interpolator = InterpolatorType::New();
resample->SetInterpolator(interpolator);
auto transform = TransformType::New();
// get transform parameters from MatrixType
TransformType::ParametersType parameters(Dimension * Dimension + Dimension);
for (unsigned int i = 0; i < Dimension; ++i)
{
for (unsigned int j = 0; j < Dimension; ++j)
{
parameters[i * Dimension + j] = matrix[i][j];
}
}
for (unsigned int i = 0; i < Dimension; ++i)
{
parameters[i + Dimension * Dimension] = matrix[i][Dimension];
}
transform->SetParameters(parameters);
resample->SetTransform(transform);
try
{
itk::WriteImage(resample->GetOutput(), outputFileName);
}
catch (const itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
itkImageFileReader.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkAffineTransform.h
itk::AffineTransform
Definition: itkAffineTransform.h:101
itk::Matrix
A templated class holding a M x N size Matrix.
Definition: itkMatrix.h:52
itkImageFileWriter.h
itk::WindowedSincInterpolateImageFunction
Use the windowed sinc function to interpolate.
Definition: itkWindowedSincInterpolateImageFunction.h:271
itkWindowedSincInterpolateImageFunction.h
itk::ResampleImageFilter
Resample an image via a coordinate transform.
Definition: itkResampleImageFilter.h:90
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itkResampleImageFilter.h
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itk::WriteImage
ITK_TEMPLATE_EXPORT void WriteImage(TImagePointer &&image, const std::string &filename, bool compress=false)
Definition: itkImageFileWriter.h:256
itk::Size::GetSize
const SizeValueType * GetSize() const
Definition: itkSize.h:169