SimpleITK  
itk::simple::ImageSeriesReader Class Reference

Read series of image files into a SimpleITK image. More...

#include <sitkImageSeriesReader.h>

+ Inheritance diagram for itk::simple::ImageSeriesReader:
+ Collaboration diagram for itk::simple::ImageSeriesReader:

Detailed Description

Read series of image files into a SimpleITK image.

For some image formats such as DICOM, images also contain associated meta-data (e.g. imaging modality, patient name etc.). By default the reader does not load this information (saves time). To load the meta-data you will need to explicitly configure the reader, MetaDataDictionaryArrayUpdateOn, and possibly specify that you also want to load the private meta-data LoadPrivateTagsOn.

Once the image series is read the meta-data is directly accessible from the reader.

Note
DICOM tags are represented as strings in the meta-data dictionary(s), therefore "0020|000D" and "0020|000d" are different when accessing the tag value. This differs from the hexadecimal numbers they represent, 0020|000D and 0020|000d are equivalent. The ITK meta-data dictionary is string based and uses lower case to represent the hexadecimal number read from disk, so 0020|000d will work as a key and 0020|000D will not be found in the dictionary (results in an exception if attempting to access). It is recommended to use lower case when setting and accessing DICOM tags.
If the pixel type for the returned image is not specified it is deduced from the first image in the series. This approach is computationally efficient and assumes that all images in a series have the same pixel type. In rare situations this is not the case, not all images have the same pixel type. If this leads to a narrowing conversion (e.g. first image pixel type is unsigned int and others are float) the returned image does not represent the data correctly. To resolve such situations, explicitly specify the expected pixel type via the SetOutputPixelType method before reading the series.
See also
itk::simple::ReadImage for the procedural interface
Examples
DicomSeriesFromArray/DicomSeriesFromArray.cs, DicomSeriesFromArray/DicomSeriesFromArray.cxx, DicomSeriesReadModifyWrite/DicomSeriesReadModifySeriesWrite.cs, DicomSeriesReadModifyWrite/DicomSeriesReadModifySeriesWrite.cxx, DicomSeriesReader/DicomSeriesReader.cs, and DicomSeriesReader/DicomSeriesReader.cxx.

Definition at line 69 of file sitkImageSeriesReader.h.

Public Types

using Self = ImageSeriesReader
 
- Public Types inherited from itk::simple::ImageReaderBase
using Self = ImageReaderBase
 
- Public Types inherited from itk::simple::ProcessObject
using Self = ProcessObject
 

Public Member Functions

Image Execute () override
 
void ForceOrthogonalDirectionOff ()
 
void ForceOrthogonalDirectionOn ()
 
const std::vector< PathType > & GetFileNames () const
 
bool GetForceOrthogonalDirection () const
 
std::string GetMetaData (unsigned int slice, const std::string &key) const
 Get the value of a meta-data dictionary entry as a string.
 
bool GetMetaDataDictionaryArrayUpdate ()
 
std::vector< std::string > GetMetaDataKeys (unsigned int slice) const
 Get the meta-data dictionary keys for a slice.
 
std::string GetName () const override
 
bool GetReverseOrder () const
 
double GetSpacingWarningRelThreshold () const
 
bool HasMetaDataKey (unsigned int slice, const std::string &key) const
 Query a meta-data dictionary for the existence of a key.
 
 ImageSeriesReader ()
 
void MetaDataDictionaryArrayUpdateOff ()
 
void MetaDataDictionaryArrayUpdateOn ()
 
void ReverseOrderOff ()
 
void ReverseOrderOn ()
 
void SetFileNames (const std::vector< PathType > &fileNames)
 
void SetForceOrthogonalDirection (bool forceOrthogonalDirection)
 
void SetMetaDataDictionaryArrayUpdate (bool metaDataDictionaryArrayUpdate)
 
void SetReverseOrder (bool reverseOrder)
 
void SetSpacingWarningRelThreshold (double spacingWarningRelThreshold)
 
std::string ToString () const override
 
 ~ImageSeriesReader () override
 
- Public Member Functions inherited from itk::simple::ImageReaderBase
virtual std::vector< std::string > GetRegisteredImageIOs () const
 Get a vector of the names of registered itk ImageIOs.
 
 ImageReaderBase ()
 
 ~ImageReaderBase () override
 
void SetOutputPixelType (PixelIDValueEnum pixelID)
 Set/Get The output PixelType of the image.
 
PixelIDValueEnum GetOutputPixelType () const
 Set/Get The output PixelType of the image.
 
virtual void SetLoadPrivateTags (bool loadPrivateTags)
 Set/Get loading private DICOM tags into Image's MetaData.
 
virtual bool GetLoadPrivateTags () const
 Set/Get loading private DICOM tags into Image's MetaData.
 
virtual void LoadPrivateTagsOn ()
 Set/Get loading private DICOM tags into Image's MetaData.
 
virtual void LoadPrivateTagsOff ()
 Set/Get loading private DICOM tags into Image's MetaData.
 
virtual void SetImageIO (const std::string &imageio)
 Set/Get name of ImageIO to use.
 
virtual std::string GetImageIO () const
 Set/Get name of ImageIO to use.
 
- Public Member Functions inherited from itk::simple::ProcessObject
virtual void Abort ()
 
virtual int AddCommand (itk::simple::EventEnum event, const std::function< void()> &func)
 Directly add a callback to observe an event.
 
virtual int AddCommand (itk::simple::EventEnum event, itk::simple::Command &cmd)
 Add a Command Object to observer the event.
 
virtual float GetProgress () const
 An Active Measurement of the progress of execution.
 
virtual bool HasCommand (itk::simple::EventEnum event) const
 Query of this object has any registered commands for event.
 
 ProcessObject ()
 
virtual void RemoveAllCommands ()
 Remove all registered commands.
 
virtual ~ProcessObject ()
 
virtual void DebugOn ()
 
virtual void DebugOff ()
 
virtual bool GetDebug () const
 
virtual void SetDebug (bool debugFlag)
 
virtual void SetNumberOfThreads (unsigned int n)
 
virtual unsigned int GetNumberOfThreads () const
 
virtual void SetNumberOfWorkUnits (unsigned int n)
 
virtual unsigned int GetNumberOfWorkUnits () const
 

Static Public Member Functions

static std::vector< PathTypeGetGDCMSeriesFileNames (const PathType &directory, const std::string &seriesID="", bool useSeriesDetails=false, bool recursive=false, bool loadSequences=false)
 Generate a sequence of filenames from a directory with a DICOM data set and a series ID.
 
static std::vector< std::string > GetGDCMSeriesIDs (const PathType &directory, bool useSeriesDetails=false)
 Get all the seriesIDs from a DICOM data set.
 
- Static Public Member Functions inherited from itk::simple::ImageReaderBase
static std::string GetImageIOFromFileName (const PathType &fileName)
 Get the automatic ImageIO from the ImageIOFactory.
 
- Static Public Member Functions inherited from itk::simple::ProcessObject
static bool GetGlobalDefaultDebug ()
 
static void GlobalDefaultDebugOff ()
 
static void GlobalDefaultDebugOn ()
 
static void SetGlobalDefaultDebug (bool debugFlag)
 
static void GlobalWarningDisplayOn ()
 
static void GlobalWarningDisplayOff ()
 
static void SetGlobalWarningDisplay (bool flag)
 
static bool GetGlobalWarningDisplay ()
 
static double GetGlobalDefaultCoordinateTolerance ()
 Access the global tolerance to determine congruent spaces.
 
static void SetGlobalDefaultCoordinateTolerance (double)
 Access the global tolerance to determine congruent spaces.
 
static double GetGlobalDefaultDirectionTolerance ()
 Access the global tolerance to determine congruent spaces.
 
static void SetGlobalDefaultDirectionTolerance (double)
 Access the global tolerance to determine congruent spaces.
 
static bool SetGlobalDefaultThreader (const std::string &threader)
 Set/Get the default threader used for process objects.
 
static std::string GetGlobalDefaultThreader ()
 Set/Get the default threader used for process objects.
 
static void SetGlobalDefaultNumberOfThreads (unsigned int n)
 
static unsigned int GetGlobalDefaultNumberOfThreads ()
 Set/Get the default threader used for process objects.
 

Protected Member Functions

template<class TImageType>
Image ExecuteInternal (itk::ImageIOBase *)
 
- Protected Member Functions inherited from itk::simple::ImageReaderBase
unsigned int GetDimensionFromImageIO (const itk::ImageIOBase *iobase, unsigned int i)
 
unsigned int GetDimensionFromImageIO (const PathType &filename, unsigned int i)
 
itk::SmartPointer< ImageIOBaseGetImageIOBase (const PathType &fileName)
 
void GetPixelIDFromImageIO (const itk::ImageIOBase *iobase, PixelIDValueType &outPixelType, unsigned int &outDimensions)
 
void GetPixelIDFromImageIO (const PathType &fileName, PixelIDValueType &outPixelType, unsigned int &outDimensions)
 
- Protected Member Functions inherited from itk::simple::ProcessObject
virtual unsigned long AddITKObserver (const itk::EventObject &, itk::Command *)
 
virtual itk::ProcessObjectGetActiveProcess ()
 
virtual void OnActiveProcessDelete ()
 
virtual void onCommandDelete (const itk::simple::Command *cmd) noexcept
 
virtual void PreUpdate (itk::ProcessObject *p)
 
virtual void RemoveITKObserver (EventCommand &e)
 
- Protected Member Functions inherited from itk::simple::NonCopyable
 NonCopyable ()=default
 
 NonCopyable (const NonCopyable &)=delete
 
NonCopyableoperator= (const NonCopyable &)=delete
 

Private Types

typedef Image(Self::* MemberFunctionType) (itk::ImageIOBase *)
 

Static Private Member Functions

static const detail::MemberFunctionFactory< MemberFunctionType > & GetMemberFunctionFactory ()
 

Private Attributes

std::vector< PathTypem_FileNames
 
std::unique_ptr< itk::ProcessObject, ProcessObjectDeleterm_Filter
 
bool m_ForceOrthogonalDirection { true }
 
bool m_MetaDataDictionaryArrayUpdate { false }
 
std::function< std::string(int, const std::string &)> m_pfGetMetaData
 
std::function< std::vector< std::string >(int)> m_pfGetMetaDataKeys
 
std::function< bool(int, const std::string &)> m_pfHasMetaDataKey
 
bool m_ReverseOrder { false }
 
double m_SpacingWarningRelThreshold { 1e-4 }
 

Friends

struct detail::MemberFunctionAddressor< MemberFunctionType >
 

Additional Inherited Members

- Static Protected Member Functions inherited from itk::simple::ProcessObject
template<class TImageType>
static TImageType::ConstPointer CastImageToITK (const Image &img)
 
template<class TPixelType, unsigned int VImageDimension, unsigned int VLength, template< typename, unsigned int > class TVector>
static Image CastITKToImage (itk::Image< TVector< TPixelType, VLength >, VImageDimension > *img)
 
template<unsigned int VImageDimension, unsigned int VLength, template< unsigned int > class TVector>
static Image CastITKToImage (itk::Image< TVector< VLength >, VImageDimension > *img)
 
template<class TImageType>
static Image CastITKToImage (TImageType *img)
 
static const itk::EventObjectGetITKEventObject (EventEnum e)
 
template<typename T>
static std::ostream & ToStringHelper (std::ostream &os, const T &v)
 
static std::ostream & ToStringHelper (std::ostream &os, const char &v)
 
static std::ostream & ToStringHelper (std::ostream &os, const signed char &v)
 
static std::ostream & ToStringHelper (std::ostream &os, const unsigned char &v)
 

Member Typedef Documentation

◆ MemberFunctionType

typedef Image(Self::* itk::simple::ImageSeriesReader::MemberFunctionType) (itk::ImageIOBase *)
private

Definition at line 258 of file sitkImageSeriesReader.h.

◆ Self

Constructor & Destructor Documentation

◆ ~ImageSeriesReader()

itk::simple::ImageSeriesReader::~ImageSeriesReader ( )
override

◆ ImageSeriesReader()

itk::simple::ImageSeriesReader::ImageSeriesReader ( )

Member Function Documentation

◆ Execute()

Image itk::simple::ImageSeriesReader::Execute ( )
overridevirtual

◆ ExecuteInternal()

template<class TImageType>
Image itk::simple::ImageSeriesReader::ExecuteInternal ( itk::ImageIOBase * )
protected

◆ ForceOrthogonalDirectionOff()

void itk::simple::ImageSeriesReader::ForceOrthogonalDirectionOff ( )

◆ ForceOrthogonalDirectionOn()

void itk::simple::ImageSeriesReader::ForceOrthogonalDirectionOn ( )

◆ GetFileNames()

const std::vector< PathType > & itk::simple::ImageSeriesReader::GetFileNames ( ) const

◆ GetForceOrthogonalDirection()

bool itk::simple::ImageSeriesReader::GetForceOrthogonalDirection ( ) const

◆ GetGDCMSeriesFileNames()

static std::vector< PathType > itk::simple::ImageSeriesReader::GetGDCMSeriesFileNames ( const PathType & directory,
const std::string & seriesID = "",
bool useSeriesDetails = false,
bool recursive = false,
bool loadSequences = false )
static

Generate a sequence of filenames from a directory with a DICOM data set and a series ID.

This method generates a sequence of filenames whose filenames point to DICOM files. The data set may contain multiple series. The seriesID string is used to select a specific series. The ordering of the filenames is based on one of several strategies, which will read all images in the directory ( assuming there is only one study/series ).

Parameters
directorySet the directory that contains the DICOM data set.
recursiveRecursively parse the input directory.
seriesIDSet the name that identifies a particular series. Default value is an empty string which will return the file names associated with the first series found in the directory. If the parameter value was obtained from a call to GDCMSeriesIDs, the value of the useSeriesDetails parameter must be the same here and in the call to GDCMSeriesIDs.
useSeriesDetailsUse additional series information such as ProtocolName and SeriesName to identify when a single SeriesUID contains multiple 3D volumes - as can occur with perfusion and DTI imaging.
loadSequencesParse any sequences in the DICOM data set. Loading DICOM files is faster when sequences are not needed.
See also
itk::GDCMSeriesFileNames
Examples
DicomSeriesFromArray/DicomSeriesFromArray.cs, DicomSeriesFromArray/DicomSeriesFromArray.cxx, DicomSeriesReadModifyWrite/DicomSeriesReadModifySeriesWrite.cs, DicomSeriesReadModifyWrite/DicomSeriesReadModifySeriesWrite.cxx, DicomSeriesReader/DicomSeriesReader.cs, and DicomSeriesReader/DicomSeriesReader.cxx.

◆ GetGDCMSeriesIDs()

static std::vector< std::string > itk::simple::ImageSeriesReader::GetGDCMSeriesIDs ( const PathType & directory,
bool useSeriesDetails = false )
static

Get all the seriesIDs from a DICOM data set.

Parameters
directoryThe directory that contains the DICOM data set
useSeriesDetailsUse additional series information such as ProtocolName and SeriesName to identify when a single SeriesUID contains multiple 3D volumes - as can occur with perfusion and DTI imaging. The parameter value must match the value used in the call to GDCMSeriesFileNames.
See also
itk::GDCMSeriesFileNames
Examples
DicomSeriesFromArray/DicomSeriesFromArray.cs, DicomSeriesFromArray/DicomSeriesFromArray.cxx, DicomSeriesReadModifyWrite/DicomSeriesReadModifySeriesWrite.cs, and DicomSeriesReadModifyWrite/DicomSeriesReadModifySeriesWrite.cxx.

◆ GetMemberFunctionFactory()

static const detail::MemberFunctionFactory< MemberFunctionType > & itk::simple::ImageSeriesReader::GetMemberFunctionFactory ( )
staticprivate

◆ GetMetaData()

std::string itk::simple::ImageSeriesReader::GetMetaData ( unsigned int slice,
const std::string & key ) const
inline

Get the value of a meta-data dictionary entry as a string.

If the key is not in the dictionary then an exception is thrown.

string types in the dictionary are returned as their native string. Other types are printed to string before returning.

Examples
DicomSeriesReadModifyWrite/DicomSeriesReadModifySeriesWrite.cxx.

Definition at line 245 of file sitkImageSeriesReader.h.

References m_pfGetMetaData.

◆ GetMetaDataDictionaryArrayUpdate()

bool itk::simple::ImageSeriesReader::GetMetaDataDictionaryArrayUpdate ( )
inline

Definition at line 101 of file sitkImageSeriesReader.h.

References m_MetaDataDictionaryArrayUpdate.

◆ GetMetaDataKeys()

std::vector< std::string > itk::simple::ImageSeriesReader::GetMetaDataKeys ( unsigned int slice) const
inline

Get the meta-data dictionary keys for a slice.

This is only valid after successful execution of this filter and when MetaDataDictionaryArrayUpdate is true. Each element in the array corresponds to a "slice" or filename read during execution.

If the slice index is out of range, an exception will be thrown.

Returns a vector of keys to the key/value entries in the file's meta-data dictionary. Iterate through with these keys to get the values.

Definition at line 223 of file sitkImageSeriesReader.h.

References m_pfGetMetaDataKeys.

◆ GetName()

std::string itk::simple::ImageSeriesReader::GetName ( ) const
inlineoverridevirtual

return user readable name of the filter

Implements itk::simple::ProcessObject.

Definition at line 84 of file sitkImageSeriesReader.h.

◆ GetReverseOrder()

bool itk::simple::ImageSeriesReader::GetReverseOrder ( ) const

◆ GetSpacingWarningRelThreshold()

double itk::simple::ImageSeriesReader::GetSpacingWarningRelThreshold ( ) const

◆ HasMetaDataKey()

bool itk::simple::ImageSeriesReader::HasMetaDataKey ( unsigned int slice,
const std::string & key ) const
inline

Query a meta-data dictionary for the existence of a key.

Examples
DicomSeriesReadModifyWrite/DicomSeriesReadModifySeriesWrite.cxx.

Definition at line 231 of file sitkImageSeriesReader.h.

References m_pfHasMetaDataKey.

◆ MetaDataDictionaryArrayUpdateOff()

void itk::simple::ImageSeriesReader::MetaDataDictionaryArrayUpdateOff ( )
inline

Definition at line 114 of file sitkImageSeriesReader.h.

References SetMetaDataDictionaryArrayUpdate().

◆ MetaDataDictionaryArrayUpdateOn()

void itk::simple::ImageSeriesReader::MetaDataDictionaryArrayUpdateOn ( )
inline

Set the value of MetaDataDictionaryArrayUpdate to true or false respectively.

Examples
DicomSeriesReadModifyWrite/DicomSeriesReadModifySeriesWrite.cs, and DicomSeriesReadModifyWrite/DicomSeriesReadModifySeriesWrite.cxx.

Definition at line 109 of file sitkImageSeriesReader.h.

References SetMetaDataDictionaryArrayUpdate().

◆ ReverseOrderOff()

void itk::simple::ImageSeriesReader::ReverseOrderOff ( )

◆ ReverseOrderOn()

void itk::simple::ImageSeriesReader::ReverseOrderOn ( )

◆ SetFileNames()

◆ SetForceOrthogonalDirection()

void itk::simple::ImageSeriesReader::SetForceOrthogonalDirection ( bool forceOrthogonalDirection)

Set/Get/OnOff whether to force orthogonal direction cosines.

Do we want to force orthogonal direction cosines? On by default. Turning it off enables proper reading of DICOM series with gantry tilt.

◆ SetMetaDataDictionaryArrayUpdate()

void itk::simple::ImageSeriesReader::SetMetaDataDictionaryArrayUpdate ( bool metaDataDictionaryArrayUpdate)
inline

Set/Get whether the meta-data dictionaries for the slices should be read. Default value is false, because of the additional computation time.

Definition at line 96 of file sitkImageSeriesReader.h.

References m_MetaDataDictionaryArrayUpdate.

Referenced by MetaDataDictionaryArrayUpdateOff(), and MetaDataDictionaryArrayUpdateOn().

◆ SetReverseOrder()

void itk::simple::ImageSeriesReader::SetReverseOrder ( bool reverseOrder)

Set/Get/OnOff whether to reverse the order of the input files.

If true, the order of the input file names will be reversed before reading. Defaults to false.

◆ SetSpacingWarningRelThreshold()

void itk::simple::ImageSeriesReader::SetSpacingWarningRelThreshold ( double spacingWarningRelThreshold)

Set the relative threshold for issuing warnings about non-uniform sampling.

◆ ToString()

std::string itk::simple::ImageSeriesReader::ToString ( ) const
overridevirtual

Print ourselves to string

Reimplemented from itk::simple::ImageReaderBase.

Friends And Related Symbol Documentation

◆ detail::MemberFunctionAddressor< MemberFunctionType >

Definition at line 258 of file sitkImageSeriesReader.h.

Member Data Documentation

◆ m_FileNames

std::vector<PathType> itk::simple::ImageSeriesReader::m_FileNames
private

Definition at line 274 of file sitkImageSeriesReader.h.

◆ m_Filter

std::unique_ptr<itk::ProcessObject, ProcessObjectDeleter> itk::simple::ImageSeriesReader::m_Filter
private

Definition at line 271 of file sitkImageSeriesReader.h.

◆ m_ForceOrthogonalDirection

bool itk::simple::ImageSeriesReader::m_ForceOrthogonalDirection { true }
private

Definition at line 280 of file sitkImageSeriesReader.h.

◆ m_MetaDataDictionaryArrayUpdate

bool itk::simple::ImageSeriesReader::m_MetaDataDictionaryArrayUpdate { false }
private

◆ m_pfGetMetaData

std::function<std::string(int, const std::string &)> itk::simple::ImageSeriesReader::m_pfGetMetaData
private

Definition at line 268 of file sitkImageSeriesReader.h.

Referenced by GetMetaData().

◆ m_pfGetMetaDataKeys

std::function<std::vector<std::string>(int)> itk::simple::ImageSeriesReader::m_pfGetMetaDataKeys
private

Definition at line 266 of file sitkImageSeriesReader.h.

Referenced by GetMetaDataKeys().

◆ m_pfHasMetaDataKey

std::function<bool(int, const std::string &)> itk::simple::ImageSeriesReader::m_pfHasMetaDataKey
private

Definition at line 267 of file sitkImageSeriesReader.h.

Referenced by HasMetaDataKey().

◆ m_ReverseOrder

bool itk::simple::ImageSeriesReader::m_ReverseOrder { false }
private

Definition at line 282 of file sitkImageSeriesReader.h.

◆ m_SpacingWarningRelThreshold

double itk::simple::ImageSeriesReader::m_SpacingWarningRelThreshold { 1e-4 }
private

Definition at line 276 of file sitkImageSeriesReader.h.


The documentation for this class was generated from the following file: