18 #include "itkImageFileReader.h" 19 #include "itkImageSeriesReader.h" 20 #include "itkImageSeriesWriter.h" 21 #include "itkCastImageFilter.h" 22 #include "itkImageFileWriter.h" 23 #include "itkImageIOBase.h" 24 #include "itkImageIOFactory.h" 25 #include "itkNiftiImageIO.h" 26 #include "itkGDCMImageIO.h" 27 #include "itkGDCMSeriesFileNames.h" 30 #include "itkNumericSeriesFileNames.h" 31 #include "itkOrientImageFilter.h" 32 #include "itkChangeInformationImageFilter.h" 34 #if ITK_VERSION_MAJOR >= 4 35 #include "gdcmUIDGenerator.h" 37 #include "gdcm/src/gdcmFile.h" 38 #include "gdcm/src/gdcmUtil.h" 44 #include "DicomIOManager.h" 46 using ImageTypeFloat3D = itk::Image< float, 3 >;
47 using TImageType = ImageTypeFloat3D;
48 using MaskType = itk::Image<unsigned int, 3>;
68 template <
class TImageType = ImageTypeFloat3D >
69 typename itk::ImageFileReader< TImageType >::Pointer
GetImageReader(
const std::string &fName,
const std::string &supportedExtensions =
".nii.gz,.nii,.dcm",
const std::string &delimitor =
",")
81 std::transform(fileExtension.begin(), fileExtension.end(), fileExtension.begin(), ::tolower);
83 auto reader = itk::ImageFileReader< TImageType >::New();
85 if (fileExtension ==
".dcm")
88 if (filesInDir.size() > 1)
90 std::cerr <<
"Trying to read DICOM file. Please use DICOMImageReader.\n";
95 if (supportedExtensions !=
"")
97 std::vector< std::string > extensions =
cbica::stringSplit(supportedExtensions, delimitor);
99 bool supportedExtensionFound =
false;
100 for (
size_t i = 0; i < extensions.size(); i++)
102 if (extensions[i] == fileExtension)
104 supportedExtensionFound =
true;
108 if (!supportedExtensionFound)
110 std::cerr <<
"Supplied file name '" << fName_wrap <<
"' doesn't have a supported extension. \nSupported Extensions: " << supportedExtensions <<
"\n";
119 if ((imageInfo.GetImageDimensions() != TImageType::ImageDimension) &&
120 !((TImageType::ImageDimension == 2) && (imageInfo.GetImageSize()[2] == 1)))
122 std::cerr <<
"Image Dimension mismatch. Return image is expected to be '" << TImageType::ImageDimension <<
123 "'D and doesn't match the image dimension read from the input file, which is '" << imageInfo.GetImageDimensions() <<
"'.\n";
127 reader->SetFileName(fName_wrap);
131 if (std::find(supportedExtsVector.begin(), supportedExtsVector.end(), fileExtension) == supportedExtsVector.end())
133 std::cerr <<
"Extension of file doesn't match the supported extensions; can't read.\n";
138 if ((fileExtension ==
".dcm") || (fileExtension ==
".dicom"))
144 else if ((fileExtension ==
".nii") || (fileExtension ==
".nii.gz"))
146 auto ioType = itk::NiftiImageIO::New();
147 ioType->SetFileName(fName_wrap);
148 reader->SetImageIO(ioType);
155 catch (itk::ExceptionObject& e)
157 std::cerr <<
"Exception caught while reading the image '" << fName_wrap <<
"': " << e.what() <<
"\n";
263 template <
class TImageType = ImageTypeFloat3D >
271 if (dirName_wrap[dirName_wrap.length() - 1] ==
'/')
272 dirName_wrap.pop_back();
281 auto dicomIO = itk::GDCMImageIO::New();
282 auto inputNames = itk::GDCMSeriesFileNames::New();
283 inputNames->SetInputDirectory(dirName_wrap);
284 inputNames->SetLoadPrivateTags(
true);
285 auto UIDs = inputNames->GetSeriesUIDs();
289 if (UIDs_unique.size() > 1)
291 std::cout <<
"Multiple DICOM series detected.\n";
294 inputNames->SetInputDirectory(dirName_wrap);
297 auto filenames = inputNames->GetInputFileNames();
299 auto seriesReader = itk::ImageSeriesReader< TImageType >::New();
300 seriesReader->SetImageIO(dicomIO);
301 seriesReader->SetFileNames(filenames);
305 seriesReader->Update();
307 catch (itk::ExceptionObject & err)
309 std::cerr <<
"Error while loading DICOM images: " << err.what() <<
"\n";
444 template <
typename ComputedImageType = ImageTypeFloat3D,
typename ExpectedImageType = ComputedImageType>
445 void WriteImage(
typename ComputedImageType::Pointer imageToWrite,
const std::string &fileName)
454 auto filter = itk::CastImageFilter<ComputedImageType, ExpectedImageType>::New();
455 filter->SetInput(imageToWrite);
458 auto writer = itk::ImageFileWriter< ExpectedImageType >::New();
461 if ((ext ==
".nii") || (ext ==
".nii.gz"))
463 writer->SetImageIO(itk::NiftiImageIO::New());
466 writer->SetInput(filter->GetOutput());
467 writer->SetFileName(fileName);
473 catch (itk::ExceptionObject &e)
475 std::cerr <<
"Error occurred while trying to write the image '" << fileName <<
"': " << e.what() <<
"\n";
502 template <
typename ComputedImageType = ImageTypeFloat3D>
503 void WriteDicomImage(
const typename ComputedImageType::Pointer imageToWrite,
const std::string &dirName)
505 auto reader =
typename itk::ImageSeriesReader< ComputedImageType >::New();
506 WriteDicomImage< ComputedImageType >(reader, imageToWrite, dirName);
528 template <
typename ComputedImageType = ImageTypeFloat3D>
529 void WriteDicomImage(
const typename itk::ImageSeriesReader< ComputedImageType >::Pointer inputImageReader,
const typename ComputedImageType::Pointer imageToWrite,
const std::string &dirName)
533 std::cout <<
"Specified directory wasn't found, creating...\n";
544 using ExpectedImageType = itk::Image< short, ComputedImageType::ImageDimension >;
545 typedef itk::CastImageFilter<ComputedImageType, ExpectedImageType> CastFilterType;
546 typename CastFilterType::Pointer castFilter = CastFilterType::New();
547 castFilter->SetInput(imageToWrite);
548 castFilter->Update();
552 auto dicomIO = itk::GDCMImageIO::New();
554 dicomIO->SetComponentType(itk::ImageIOBase::IOComponentType::SHORT);
556 auto seriesWriter = itk::ImageSeriesWriter< ExpectedImageType, itk::Image<typename ExpectedImageType::PixelType, 2> >::New();
558 auto namesGenerator = itk::NumericSeriesFileNames::New();
560 auto start = imageToWrite->GetLargestPossibleRegion().GetIndex();
561 auto size = imageToWrite->GetLargestPossibleRegion().GetSize();
562 namesGenerator->SetSeriesFormat((dirName +
"/image%03d.dcm").c_str());
563 namesGenerator->SetStartIndex(start[2]);
564 namesGenerator->SetEndIndex(start[2] + size[2] - 1);
565 namesGenerator->SetIncrementIndex(1);
567 seriesWriter->SetInput(castFilter->GetOutput());
568 seriesWriter->SetImageIO(dicomIO);
569 seriesWriter->SetFileNames(namesGenerator->GetFileNames());
571 typename itk::ImageSeriesReader< ComputedImageType >::DictionaryArrayType outputArray;
572 if (inputImageReader.IsNull() || (inputImageReader->GetImageIO() == NULL))
584 typename ExpectedImageType::IndexType index;
587 for (
size_t i = 0; i < imageToWrite->GetLargestPossibleRegion().GetSize()[2]; i++)
589 auto dict =
new typename itk::ImageSeriesReader< ComputedImageType >::DictionaryType;
590 typename ExpectedImageType::PointType position;
592 imageToWrite->TransformIndexToPhysicalPoint(index, position);
593 itk::EncapsulateMetaData<std::string>(*dict,
"0020|0032", std::to_string(position[0]) +
"\\" + std::to_string(position[1]) +
"\\" + std::to_string(position[2]));
594 itk::EncapsulateMetaData<std::string>(*dict,
"0020|1041", std::to_string(position[0]) +
"\\" + std::to_string(position[1]) +
"\\" + std::to_string(position[2]));
613 itk::EncapsulateMetaData<std::string>(*dict,
"0018|0050", std::to_string(imageToWrite->GetSpacing()[2]));
614 itk::EncapsulateMetaData<std::string>(*dict,
"0018|0088", std::to_string(imageToWrite->GetSpacing()[2]));
615 itk::EncapsulateMetaData<std::string>(*dict,
"0028|0030", std::to_string(imageToWrite->GetSpacing()[0]) +
"\\" + std::to_string(imageToWrite->GetSpacing()[1]));
621 outputArray.push_back(dict);
624 seriesWriter->SetMetaDataDictionaryArray(&outputArray);
628 dicomIO->SetMetaDataDictionary(inputImageReader->GetMetaDataDictionary());
629 seriesWriter->SetMetaDataDictionaryArray(inputImageReader->GetMetaDataDictionaryArray());
634 seriesWriter->Write();
636 catch (itk::ExceptionObject &e)
638 std::cerr <<
"Error occurred while trying to write the image '" << dirName <<
"': " << e.what() <<
"\n";
661 template <
class TImageType = ImageTypeFloat3D >
662 typename TImageType::Pointer
ReadImage(
const std::string &fName,
const std::string &supportedExtensions =
".nii.gz,.nii,.dcm",
const std::string &delimitor =
",")
664 if (cbica::IsDicom(fName))
668 bool loadstatus = dcmSeriesReader.
LoadDicom();
678 return GetImageReader< TImageType >(fName, supportedExtensions, delimitor)->GetOutput();
688 template<
class TImageType >
691 auto orienter = itk::OrientImageFilter<TImageType, TImageType>::New();
692 orienter->UseImageDirectionOn();
693 orienter->SetDesiredCoordinateOrientation(itk::SpatialOrientation::ITK_COORDINATE_ORIENTATION_RAI);
694 orienter->SetInput(inputImage);
697 return orienter->GetOutput();
716 template<
class TImageType >
717 typename TImageType::Pointer
ReadImageWithOrientFix(
const std::string &fName,
const std::string &supportedExtensions =
".nii.gz,.nii",
const std::string &delimitor =
",")
720 std::transform(extension.begin(), extension.end(), extension.begin(), ::tolower);
727 return GetImageWithOrientFix<TImageType>(GetImageReader< TImageType >(fName, supportedExtensions, delimitor)->GetOutput());
747 template <
class TImageType = ImageTypeFloat3D >
748 typename TImageType::Pointer
GetImage(
const std::string &fName,
const std::string &supportedExtensions =
".nii.gz,.nii",
const std::string &delimitor =
",")
750 return ReadImage< TImageType >(fName, supportedExtensions, delimitor);
bool isDir(const std::string &path)
Return True if path is an existing directory.
void SetDirectoryPath(std::string path)
set the input directory containing dicom series
Definition: DicomIOManager.hxx:18
TImageType::Pointer ReadImageWithOrientFix(const std::string &fName, const std::string &supportedExtensions=".nii.gz,.nii", const std::string &delimitor=",")
The reads the image according to the appropriate extension and outputs the result in ITK's RAI orient...
Definition: cbicaITKSafeImageIO.h:717
void WriteImage(typename ComputedImageType::Pointer imageToWrite, const std::string &fileName)
Get the itk::Image from input dir name.
Definition: cbicaITKSafeImageIO.h:445
std::vector< std::string > stringSplit(const std::string &str, const std::string &delim)
Splits the string.
TImageType::Pointer GetImageWithOrientFix(const typename TImageType::Pointer inputImage)
This is an inline function used to correct the orientation for correct visualization.
Definition: cbicaITKSafeImageIO.h:689
Definition: DicomIOManager.h:24
T::Pointer GetITKImage()
get the read dicom data as 3D float ITK image
Definition: DicomIOManager.hxx:308
Reads any image from file name and generates relevant data.
Definition: cbicaITKImageInfo.h:36
TImageType::Pointer ReadImage(const std::string &fName, const std::string &supportedExtensions=".nii.gz,.nii,.dcm", const std::string &delimitor=",")
Get the itk::Image from input file name.
Definition: cbicaITKSafeImageIO.h:662
Some basic utility functions.
bool LoadDicom()
load dicom data
Definition: DicomIOManager.hxx:23
std::vector< std::string > filesInDirectory(const std::string &dirName, bool returnFullPath=true)
Find all files inside a directory.
std::string normPath(const std::string &path)
Normalize a pathname by collapsing redundant separators and up-level references.
std::string getFilenameExtension(const std::string &filename, bool checkFile=true)
Gets the extension of the supplied file name using splitFileName()
Some basic utility functions.
itk::ImageSeriesReader< TImageType >::Pointer GetDicomImageReader(const std::string &dirName)
Returns the unique series IDs in the specified directory.
Definition: cbicaITKSafeImageIO.h:264
itk::ImageFileReader< TImageType >::Pointer GetImageReader(const std::string &fName, const std::string &supportedExtensions=".nii.gz,.nii,.dcm", const std::string &delimitor=",")
Get the itk::ImageFileReader from input file name.
Definition: cbicaITKSafeImageIO.h:69
std::string getFilenamePath(const std::string &filename, bool checkFile=true)
Gets the path of the supplied file name using splitFileName()
Declaration of the ImageInfo class.
TImageType::Pointer GetImage(const std::string &fName, const std::string &supportedExtensions=".nii.gz,.nii", const std::string &delimitor=",")
Get the itk::Image from input file name.
Definition: cbicaITKSafeImageIO.h:748
bool createDir(const std::string &dir_name)
Create a directory.
std::vector< TDataType > GetUniqueElements(const std::vector< TDataType > &inputVector)
Find the unique elements in a vector.
Definition: cbicaUtilities.h:711