SimpleITK  2.1.0
sitkImageSeriesReader.h
Go to the documentation of this file.
1 /*=========================================================================
2 *
3 * Copyright NumFOCUS
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License");
6 * you may not use this file except in compliance with the License.
7 * You may obtain a copy of the License at
8 *
9 * http://www.apache.org/licenses/LICENSE-2.0.txt
10 *
11 * Unless required by applicable law or agreed to in writing, software
12 * distributed under the License is distributed on an "AS IS" BASIS,
13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 * See the License for the specific language governing permissions and
15 * limitations under the License.
16 *
17 *=========================================================================*/
18 #ifndef sitkImageSeriesReader_h
19 #define sitkImageSeriesReader_h
20 
21 #include "sitkMacro.h"
22 #include "sitkImage.h"
23 #include "sitkImageReaderBase.h"
25 
26 namespace itk {
27  namespace simple {
28 
46  : public ImageReaderBase
47  {
48  public:
50 
51  ~ImageSeriesReader() override;
52 
54 
56  std::string ToString() const override;
57 
59  std::string GetName() const override { return std::string("ImageSeriesReader"); }
60 
61 
67  SITK_RETURN_SELF_TYPE_HEADER SetMetaDataDictionaryArrayUpdate ( bool metaDataDictionaryArrayUpdate )
68  { this->m_MetaDataDictionaryArrayUpdate = metaDataDictionaryArrayUpdate; return *this; }
69  bool GetMetaDataDictionaryArrayUpdate() { return this->m_MetaDataDictionaryArrayUpdate; }
70 
71 
73  SITK_RETURN_SELF_TYPE_HEADER MetaDataDictionaryArrayUpdateOn() { return this->SetMetaDataDictionaryArrayUpdate(true); }
74  SITK_RETURN_SELF_TYPE_HEADER MetaDataDictionaryArrayUpdateOff() { return this->SetMetaDataDictionaryArrayUpdate(false); }
75 
76 
100  static std::vector<std::string> GetGDCMSeriesFileNames( const std::string &directory,
101  const std::string &seriesID = "",
102  bool useSeriesDetails = false,
103  bool recursive = false,
104  bool loadSequences = false );
105 
111  static std::vector<std::string> GetGDCMSeriesIDs( const std::string &directory );
112 
113  SITK_RETURN_SELF_TYPE_HEADER SetFileNames ( const std::vector<std::string> &fileNames );
114  const std::vector<std::string> &GetFileNames() const;
115 
116  Image Execute() override;
117 
132  std::vector<std::string> GetMetaDataKeys( unsigned int slice ) const { return this->m_pfGetMetaDataKeys(slice); }
133 
136  bool HasMetaDataKey( unsigned int slice, const std::string &key ) const { return this->m_pfHasMetaDataKey(slice, key); }
137 
146  std::string GetMetaData( unsigned int slice, const std::string &key ) const { return this->m_pfGetMetaData(slice, key); }
147 
148 
149  protected:
150 
151  template <class TImageType> Image ExecuteInternal ( itk::ImageIOBase * );
152 
153  private:
154 
155  // function pointer type
156  typedef Image (Self::*MemberFunctionType)( itk::ImageIOBase * );
157 
158  // friend to get access to executeInternal member
159  friend struct detail::MemberFunctionAddressor<MemberFunctionType>;
160  std::unique_ptr<detail::MemberFunctionFactory<MemberFunctionType> > m_MemberFactory;
161 
162 
163  std::function<std::vector<std::string>(int)> m_pfGetMetaDataKeys;
164  std::function<bool(int, const std::string &)> m_pfHasMetaDataKey;
165  std::function<std::string(int, const std::string &)> m_pfGetMetaData;
166 
167  // Holder of process object for active measurements
169 
170 
171  std::vector<std::string> m_FileNames;
172 
174  };
175 
193  SITKIO_EXPORT Image ReadImage ( const std::vector<std::string> &fileNames,
194  PixelIDValueEnum outputPixelType=sitkUnknown,
195  const std::string &imageIO="");
196  }
197 }
198 
199 #endif
itk::simple::Image
The Image class for SimpleITK.
Definition: sitkImage.h:75
itk::ImageIOBase
itk::simple::ImageSeriesReader::m_pfHasMetaDataKey
std::function< bool(int, const std::string &)> m_pfHasMetaDataKey
Definition: sitkImageSeriesReader.h:164
itk::simple::ImageSeriesReader::GetMetaDataDictionaryArrayUpdate
bool GetMetaDataDictionaryArrayUpdate()
Definition: sitkImageSeriesReader.h:69
itk::simple::ImageSeriesReader::m_pfGetMetaData
std::function< std::string(int, const std::string &)> m_pfGetMetaData
Definition: sitkImageSeriesReader.h:165
itk::simple::detail::MemberFunctionAddressor
Definition: sitkDetail.h:32
itk::simple::ImageSeriesReader::GetMetaDataKeys
std::vector< std::string > GetMetaDataKeys(unsigned int slice) const
Get the meta-data dictionary keys for a slice.
Definition: sitkImageSeriesReader.h:132
itk::simple::ImageSeriesReader::GetMetaData
std::string GetMetaData(unsigned int slice, const std::string &key) const
Get the value of a meta-data dictionary entry as a string.
Definition: sitkImageSeriesReader.h:146
itk::simple::ImageSeriesReader::HasMetaDataKey
bool HasMetaDataKey(unsigned int slice, const std::string &key) const
Query a meta-data dictionary for the existence of a key.
Definition: sitkImageSeriesReader.h:136
SITKIO_EXPORT
#define SITKIO_EXPORT
Definition: sitkIO.h:33
sitkImage.h
itk::simple::PixelIDValueEnum
PixelIDValueEnum
Enumerated values of pixelIDs.
Definition: sitkPixelIDValues.h:90
sitkMemberFunctionFactory.h
itk::simple::ImageSeriesReader::m_pfGetMetaDataKeys
std::function< std::vector< std::string >int)> m_pfGetMetaDataKeys
Definition: sitkImageSeriesReader.h:163
sitkMacro.h
sitkImageReaderBase.h
itk::simple::ImageSeriesReader
Read series of image files into a SimpleITK image.
Definition: sitkImageSeriesReader.h:45
itk::simple::ImageSeriesReader::MetaDataDictionaryArrayUpdateOff
Self & MetaDataDictionaryArrayUpdateOff()
Definition: sitkImageSeriesReader.h:74
itk::simple::ReadImage
SITKIO_EXPORT Image ReadImage(const std::string &filename, PixelIDValueEnum outputPixelType=sitkUnknown, const std::string &imageIO="")
ReadImage is a procedural interface to the ImageFileReader class which is convenient for most image r...
itk::simple::ImageSeriesReader::SetMetaDataDictionaryArrayUpdate
Self & SetMetaDataDictionaryArrayUpdate(bool metaDataDictionaryArrayUpdate)
Definition: sitkImageSeriesReader.h:67
itk::simple::ImageReaderBase
An abstract base class for image readers.
Definition: sitkImageReaderBase.h:39
itk::simple::ImageSeriesReader::m_MetaDataDictionaryArrayUpdate
bool m_MetaDataDictionaryArrayUpdate
Definition: sitkImageSeriesReader.h:173
itk::simple::ImageSeriesReader::m_MemberFactory
std::unique_ptr< detail::MemberFunctionFactory< MemberFunctionType > > m_MemberFactory
Definition: sitkImageSeriesReader.h:160
itk::simple::ImageSeriesReader::MetaDataDictionaryArrayUpdateOn
Self & MetaDataDictionaryArrayUpdateOn()
Definition: sitkImageSeriesReader.h:73
itk
itk::simple::ImageSeriesReader::m_Filter
itk::ProcessObject * m_Filter
Definition: sitkImageSeriesReader.h:168
itk::simple::ImageSeriesReader::m_FileNames
std::vector< std::string > m_FileNames
Definition: sitkImageSeriesReader.h:171
itk::ProcessObject
itk::simple::ProcessObject
Base class for SimpleITK classes based on ProcessObject.
Definition: sitkProcessObject.h:51
itk::simple::ImageSeriesReader::GetName
std::string GetName() const override
Definition: sitkImageSeriesReader.h:59
itk::simple::sitkUnknown
@ sitkUnknown
Definition: sitkPixelIDValues.h:91