SimpleITK  
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Types | Private Attributes | Friends | List of all members
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
DicomSeriesReader/DicomSeriesReader.cxx.

Definition at line 68 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
 Set/Get The output PixelType of the image. More...
 
const std::vector< PathType > & GetFileNames () const
 
std::string GetMetaData (unsigned int slice, const std::string &key) const
 Get the value of a meta-data dictionary entry as a string. More...
 
bool GetMetaDataDictionaryArrayUpdate ()
 
std::vector< std::string > GetMetaDataKeys (unsigned int slice) const
 Get the meta-data dictionary keys for a slice. More...
 
std::string GetName () const override
 
double GetSpacingWarningRelThreshold () const
 
bool HasMetaDataKey (unsigned int slice, const std::string &key) const
 Query a meta-data dictionary for the existence of a key. More...
 
 ImageSeriesReader ()
 
SelfMetaDataDictionaryArrayUpdateOff ()
 
SelfMetaDataDictionaryArrayUpdateOn ()
 
SelfSetFileNames (const std::vector< PathType > &fileNames)
 
SelfSetMetaDataDictionaryArrayUpdate (bool metaDataDictionaryArrayUpdate)
 
SelfSetSpacingWarningRelThreshold (double spacingWarningRelThreshold)
 
std::string ToString () const override
 
 ~ImageSeriesReader () override
 
- Public Member Functions inherited from itk::simple::ImageReaderBase
 ImageReaderBase ()
 
 ~ImageReaderBase () override
 
SelfSetOutputPixelType (PixelIDValueEnum pixelID)
 Set/Get The output PixelType of the image. More...
 
PixelIDValueEnum GetOutputPixelType () const
 Set/Get The output PixelType of the image. More...
 
virtual std::vector< std::string > GetRegisteredImageIOs () const
 Get a vector of the names of registered itk ImageIOs. More...
 
virtual SelfSetLoadPrivateTags (bool loadPrivateTags)
 Set/Get loading private DICOM tags into Image's MetaData. More...
 
virtual bool GetLoadPrivateTags () const
 Set/Get The output PixelType of the image. More...
 
virtual void LoadPrivateTagsOn ()
 Set/Get The output PixelType of the image. More...
 
virtual void LoadPrivateTagsOff ()
 Set/Get The output PixelType of the image. More...
 
virtual SelfSetImageIO (const std::string &imageio)
 Set/Get name of ImageIO to use. More...
 
virtual std::string GetImageIO () const
 Set/Get The output PixelType of the image. More...
 
- 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. More...
 
virtual int AddCommand (itk::simple::EventEnum event, itk::simple::Command &cmd)
 Add a Command Object to observer the event. More...
 
virtual float GetProgress () const
 An Active Measurement of the progress of execution. More...
 
virtual bool HasCommand (itk::simple::EventEnum event) const
 Query of this object has any registered commands for event. More...
 
 ProcessObject ()
 
virtual void RemoveAllCommands ()
 Remove all registered commands. More...
 
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. More...
 
static std::vector< std::string > GetGDCMSeriesIDs (const PathType &directory, bool useSeriesDetails=false)
 Get all the seriesIDs from a DICOM data set. More...
 
- Static Public Member Functions inherited from itk::simple::ImageReaderBase
static std::string GetImageIOFromFileName (const PathType &fileName)
 Get the automatic ImageIO from the ImageIOFactory. More...
 
- 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. More...
 
static void SetGlobalDefaultCoordinateTolerance (double)
 Access the global tolerance to determine congruent spaces. More...
 
static double GetGlobalDefaultDirectionTolerance ()
 Access the global tolerance to determine congruent spaces. More...
 
static void SetGlobalDefaultDirectionTolerance (double)
 Access the global tolerance to determine congruent spaces. More...
 
static bool SetGlobalDefaultThreader (const std::string &threader)
 Set/Get the default threader used for process objects. More...
 
static std::string GetGlobalDefaultThreader ()
 Set/Get the default threader used for process objects. More...
 
static void SetGlobalDefaultNumberOfThreads (unsigned int n)
 
static unsigned int GetGlobalDefaultNumberOfThreads ()
 Set/Get the default threader used for process objects. More...
 

Protected Member Functions

template<class TImageType >
Image ExecuteInternal (itk::ImageIOBase *)
 
- Protected Member Functions inherited from itk::simple::ImageReaderBase
itk::SmartPointer< ImageIOBaseGetImageIOBase (const PathType &fileName)
 Set/Get The output PixelType of the image. More...
 
void GetPixelIDFromImageIO (const PathType &fileName, PixelIDValueType &outPixelType, unsigned int &outDimensions)
 Set/Get The output PixelType of the image. More...
 
void GetPixelIDFromImageIO (const itk::ImageIOBase *iobase, PixelIDValueType &outPixelType, unsigned int &outDimensions)
 Set/Get The output PixelType of the image. More...
 
unsigned int GetDimensionFromImageIO (const PathType &filename, unsigned int i)
 Set/Get The output PixelType of the image. More...
 
unsigned int GetDimensionFromImageIO (const itk::ImageIOBase *iobase, unsigned int i)
 Set/Get The output PixelType of the image. More...
 
- 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 *)
 

Private Attributes

std::vector< PathTypem_FileNames
 
itk::ProcessObjectm_Filter {nullptr}
 
std::unique_ptr< detail::MemberFunctionFactory< MemberFunctionType > > m_MemberFactory
 
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
 
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 231 of file sitkImageSeriesReader.h.

◆ Self

Definition at line 71 of file sitkImageSeriesReader.h.

Constructor & Destructor Documentation

◆ ~ImageSeriesReader()

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

◆ ImageSeriesReader()

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

Member Function Documentation

◆ Execute()

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

Set/Get The output PixelType of the image.

By default the value is sitkUnknown, which enable the output pixel type to be same as the file. If the pixel type is specified then the itk::ConvertPixelBuffer will be used to convert the pixels.

Implements itk::simple::ImageReaderBase.

Examples
DicomSeriesReader/DicomSeriesReader.cxx.

◆ ExecuteInternal()

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

◆ GetFileNames()

const std::vector<PathType>& itk::simple::ImageSeriesReader::GetFileNames ( ) 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
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

◆ 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.

Definition at line 218 of file sitkImageSeriesReader.h.

◆ GetMetaDataDictionaryArrayUpdate()

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

Definition at line 101 of file sitkImageSeriesReader.h.

◆ 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 196 of file sitkImageSeriesReader.h.

◆ GetName()

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

return user readable name of the filter

Implements itk::simple::ProcessObject.

Definition at line 83 of file sitkImageSeriesReader.h.

◆ 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.

Definition at line 204 of file sitkImageSeriesReader.h.

◆ MetaDataDictionaryArrayUpdateOff()

Self& itk::simple::ImageSeriesReader::MetaDataDictionaryArrayUpdateOff ( )
inline

Definition at line 114 of file sitkImageSeriesReader.h.

◆ MetaDataDictionaryArrayUpdateOn()

Self& itk::simple::ImageSeriesReader::MetaDataDictionaryArrayUpdateOn ( )
inline

Set the value of MetaDataDictionaryArrayUpdate to true or false respectively.

Definition at line 109 of file sitkImageSeriesReader.h.

◆ SetFileNames()

Self& itk::simple::ImageSeriesReader::SetFileNames ( const std::vector< PathType > &  fileNames)

◆ SetMetaDataDictionaryArrayUpdate()

Self& 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 95 of file sitkImageSeriesReader.h.

◆ SetSpacingWarningRelThreshold()

Self& 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 Function Documentation

◆ detail::MemberFunctionAddressor< MemberFunctionType >

Definition at line 234 of file sitkImageSeriesReader.h.

Member Data Documentation

◆ m_FileNames

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

Definition at line 246 of file sitkImageSeriesReader.h.

◆ m_Filter

itk::ProcessObject* itk::simple::ImageSeriesReader::m_Filter {nullptr}
private

Definition at line 243 of file sitkImageSeriesReader.h.

◆ m_MemberFactory

std::unique_ptr<detail::MemberFunctionFactory<MemberFunctionType> > itk::simple::ImageSeriesReader::m_MemberFactory
private

Definition at line 235 of file sitkImageSeriesReader.h.

◆ m_MetaDataDictionaryArrayUpdate

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

Definition at line 250 of file sitkImageSeriesReader.h.

◆ m_pfGetMetaData

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

Definition at line 240 of file sitkImageSeriesReader.h.

◆ m_pfGetMetaDataKeys

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

Definition at line 238 of file sitkImageSeriesReader.h.

◆ m_pfHasMetaDataKey

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

Definition at line 239 of file sitkImageSeriesReader.h.

◆ m_SpacingWarningRelThreshold

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

Definition at line 248 of file sitkImageSeriesReader.h.


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