ITK  6.0.0
Insight Toolkit
Examples/IO/DicomSeriesReadImageWrite.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.
*
*=========================================================================*/
//====================================================================
//
//
// WARNING WARNING WARNING WARNING WARNING WARNING !!!!
// WARNING WARNING WARNING WARNING WARNING WARNING !!!!
//
//
//
// The DICOM classes illustrated in this example are OBSOLETE and are
// scheduled for being removed from the toolkit.
//
// Please refer to DicomSeriesReadImageWrite2.cxx instead, where the GDCM
// classes are used. GDCM classes are the ones currently recommended for
// performing reading and writing of DICOM images.
//
//
//
// WARNING WARNING WARNING WARNING WARNING WARNING !!!!
// WARNING WARNING WARNING WARNING WARNING WARNING !!!!
//
//
//====================================================================
// Software Guide : BeginLatex
//
// Probably the most common representation of datasets in clinical
// applications is the one that uses sets of DICOM slices in order to compose
// tridimensional images. This is the case for CT, MRI and PET scanners. It
// is very common therefore for image analysts to have to process
// volumetric images that are stored in the form of a set of DICOM files
// belonging to a common DICOM series.
//
// The following example illustrates how to use ITK functionalities in order
// to read a DICOM series into a volume and then save this volume in another
// file format.
//
// Software Guide : EndLatex
#include "itkDICOMImageIO2.h"
#include "itkDICOMSeriesFileNames.h"
int
main(int argc, char * argv[])
{
if (argc < 3)
{
std::cerr << "Usage: " << argv[0]
<< " DicomDirectory outputFileName [seriesName]" << std::endl;
return EXIT_FAILURE;
}
using ImageType = itk::Image<short, 3>;
// Get the DICOM filenames from the directory
nameGenerator->SetDirectory(argv[1]);
try
{
using seriesIdContainer = std::vector<std::string>;
const seriesIdContainer & seriesUID = nameGenerator->GetSeriesUIDs();
seriesIdContainer::const_iterator seriesItr = seriesUID.begin();
seriesIdContainer::const_iterator seriesEnd = seriesUID.end();
std::cout << std::endl << "The directory: " << std::endl;
std::cout << std::endl << argv[1] << std::endl << std::endl;
std::cout << "Contains the following DICOM Series: ";
std::cout << std::endl << std::endl;
while (seriesItr != seriesEnd)
{
std::cout << seriesItr->c_str() << std::endl;
++seriesItr;
}
std::cout << std::endl << std::endl;
std::cout << "Now reading series: " << std::endl << std::endl;
using fileNamesContainer = std::vector<std::string>;
fileNamesContainer fileNames;
if (argc < 4) // If no optional third argument
{
std::cout << seriesUID.begin()->c_str() << std::endl;
fileNames = nameGenerator->GetFileNames();
}
else
{
std::cout << argv[3] << std::endl;
fileNames = nameGenerator->GetFileNames(argv[3]);
}
std::cout << std::endl << std::endl;
auto reader = ReaderType::New();
reader->SetFileNames(fileNames);
reader->SetImageIO(dicomIO);
try
{
reader->Update();
}
catch (const itk::ExceptionObject & ex)
{
std::cout << ex << std::endl;
return EXIT_FAILURE;
}
using WriterType = itk::ImageFileWriter<ImageType>;
auto writer = WriterType::New();
std::cout << "Writing the image as " << std::endl << std::endl;
std::cout << argv[2] << std::endl << std::endl;
writer->SetFileName(argv[2]);
writer->SetInput(reader->GetOutput());
try
{
writer->Update();
}
catch (const itk::ExceptionObject & ex)
{
std::cout << ex;
return EXIT_FAILURE;
}
}
catch (const itk::ExceptionObject & ex)
{
std::cout << ex;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::ImageSeriesReader
Data source that reads image data from a series of disk files.
Definition: itkImageSeriesReader.h:45
itkImageSeriesReader.h
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:90
itkImageFileWriter.h
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()