ITK  6.0.0
Insight Toolkit
SphinxExamples/src/IO/ImageBase/GenerateSlicesFromVolume/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 "itkImage.h"
int
main(int argc, char * argv[])
{
if (argc != 3)
{
std::cerr << "Usage: " << std::endl;
std::cerr << argv[0];
std::cerr << " <InputFileName> <OutputFileName>";
std::cerr << std::endl;
return EXIT_FAILURE;
}
const char * inputFileName = argv[1];
const char * outputFileName = argv[2];
std::string format = std::string(outputFileName) + std::string("-%d.png");
constexpr unsigned int Dimension = 3;
using PixelType = unsigned char;
using InputImageType = itk::Image<PixelType, Dimension>;
const auto input = itk::ReadImage<InputImageType>(inputFileName);
using OutputPixelType = unsigned char;
using RescaleImageType = itk::Image<OutputPixelType, Dimension>;
auto rescale = RescaleFilterType::New();
rescale->SetInput(input);
rescale->SetOutputMinimum(0);
rescale->SetOutputMaximum(255);
rescale->UpdateLargestPossibleRegion();
InputImageType::RegionType region = input->GetLargestPossibleRegion();
fnames->SetStartIndex(0);
fnames->SetEndIndex(size[2] - 1);
fnames->SetIncrementIndex(1);
fnames->SetSeriesFormat(format.c_str());
using OutputImageType = itk::Image<OutputPixelType, 2>;
auto writer = WriterType::New();
writer->SetInput(rescale->GetOutput());
writer->SetFileNames(fnames->GetFileNames());
try
{
writer->Update();
}
catch (const itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
itkImageSeriesWriter.h
itk::ImageSeriesWriter
Writes image data to a series of data files.
Definition: itkImageSeriesWriter.h:85
itkImageFileReader.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itk::SmartPointer< Self >
itkNumericSeriesFileNames.h
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkRescaleIntensityImageFilter.h
itk::NumericSeriesFileNames::New
static Pointer New()
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
itk::Size::GetSize
const SizeValueType * GetSize() const
Definition: itkSize.h:169