ITK  5.2.0
Insight Toolkit
itkImageSeriesReader.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 itkImageSeriesReader_h
19 #define itkImageSeriesReader_h
20 #include "ITKIOImageBaseExport.h"
21 
22 #include "itkSize.h"
23 #include <vector>
24 #include <string>
25 #include "itkMetaDataDictionary.h"
26 #include "itkImageFileReader.h"
27 
28 namespace itk
29 {
44 template <typename TOutputImage>
45 class ITK_TEMPLATE_EXPORT ImageSeriesReader : public ImageSource<TOutputImage>
46 {
47 public:
48  ITK_DISALLOW_COPY_AND_MOVE(ImageSeriesReader);
49 
54 
56  itkNewMacro(Self);
57 
59  itkTypeMacro(ImageSeriesReader, ImageSource);
60 
62  using SizeType = typename TOutputImage::SizeType;
63 
66 
69 
71  using OutputImagePixelType = typename TOutputImage::PixelType;
72 
76  using DictionaryArrayType = std::vector<DictionaryRawPointer>;
78 
79  using FileNamesContainer = std::vector<std::string>;
80 
83  void
85  {
86  if (m_FileNames != name)
87  {
88  m_FileNames = name;
89  this->Modified();
90  }
91  }
93 
94  const FileNamesContainer &
95  GetFileNames() const
96  {
97  return m_FileNames;
98  }
99 
102  void
103  SetFileName(std::string const & name)
104  {
105  m_FileNames.clear();
106  m_FileNames.push_back(name);
107  this->Modified();
108  }
110 
112  void
113  AddFileName(std::string const & name)
114  {
115  m_FileNames.push_back(name);
116  this->Modified();
117  }
119 
122  itkSetMacro(ReverseOrder, bool);
123  itkGetConstMacro(ReverseOrder, bool);
124  itkBooleanMacro(ReverseOrder);
126 
129  itkSetMacro(ForceOrthogonalDirection, bool);
130  itkGetConstMacro(ForceOrthogonalDirection, bool);
131  itkBooleanMacro(ForceOrthogonalDirection);
133 
138  itkSetObjectMacro(ImageIO, ImageIOBase);
139  itkGetModifiableObjectMacro(ImageIO, ImageIOBase);
141 
151  itkSetMacro(MetaDataDictionaryArrayUpdate, bool);
152  itkGetConstMacro(MetaDataDictionaryArrayUpdate, bool);
153  itkBooleanMacro(MetaDataDictionaryArrayUpdate);
154 
157  void
158  GenerateOutputInformation() override;
159 
165  void
166  EnlargeOutputRequestedRegion(DataObject * output) override;
167 
170  DictionaryArrayRawPointer
171  GetMetaDataDictionaryArray() const;
172 
174  itkSetMacro(UseStreaming, bool);
175  itkGetConstReferenceMacro(UseStreaming, bool);
176  itkBooleanMacro(UseStreaming);
178 
180  itkSetMacro(SpacingWarningRelThreshold, double);
181  itkGetConstMacro(SpacingWarningRelThreshold, double);
183 
184 protected:
186  : m_ImageIO(nullptr)
187 
188  {}
189  ~ImageSeriesReader() override;
190  void
191  PrintSelf(std::ostream & os, Indent indent) const override;
192 
194  void
195  GenerateData() override;
196 
199 
201  bool m_ReverseOrder{ false };
202 
204  bool m_ForceOrthogonalDirection{ true };
205 
208 
214  unsigned int m_NumberOfDimensionsInImage{ 0 };
215 
219 
220  bool m_UseStreaming{ true };
221 
222  bool m_SpacingDefined{ false };
223 
224  double m_SpacingWarningRelThreshold{ 1e-4 };
225 
226 private:
228 
229  int
230  ComputeMovingDimensionIndex(ReaderType * reader);
231 
234 
236  bool m_MetaDataDictionaryArrayUpdate{ true };
237 };
238 } // namespace itk
239 
240 #ifndef ITK_MANUAL_INSTANTIATION
241 # include "itkImageSeriesReader.hxx"
242 #endif
243 
244 #endif // itkImageSeriesReader_h
itk::ImageSeriesReader::SetFileNames
void SetFileNames(const FileNamesContainer &name)
Definition: itkImageSeriesReader.h:84
itk::ImageIOBase
Abstract superclass defines image IO interface.
Definition: itkImageIOBase.h:77
itk::ImageSeriesReader::ImageSeriesReader
ImageSeriesReader()
Definition: itkImageSeriesReader.h:185
itk::ImageSeriesReader::m_MetaDataDictionaryArray
DictionaryArrayType m_MetaDataDictionaryArray
Definition: itkImageSeriesReader.h:218
itk::ImageSeriesReader::SetFileName
void SetFileName(std::string const &name)
Definition: itkImageSeriesReader.h:103
itkImageFileReader.h
itk::ImageSeriesReader
Data source that reads image data from a series of disk files.
Definition: itkImageSeriesReader.h:45
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::ImageSeriesReader::GetFileNames
const FileNamesContainer & GetFileNames() const
Definition: itkImageSeriesReader.h:95
itk::ImageSeriesReader::m_MetaDataDictionaryArrayMTime
TimeStamp m_MetaDataDictionaryArrayMTime
Definition: itkImageSeriesReader.h:233
itk::SmartPointer< Self >
itk::ImageSeriesReader::SizeType
typename TOutputImage::SizeType SizeType
Definition: itkImageSeriesReader.h:62
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::ImageSeriesReader::m_ImageIO
ImageIOBase::Pointer m_ImageIO
Definition: itkImageSeriesReader.h:198
itk::ImageSeriesReader::ImageRegionType
typename TOutputImage::RegionType ImageRegionType
Definition: itkImageSeriesReader.h:68
itk::MetaDataDictionary
Provides a mechanism for storing a collection of arbitrary data types.
Definition: itkMetaDataDictionary.h:53
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ImageSeriesReader::IndexType
typename TOutputImage::IndexType IndexType
Definition: itkImageSeriesReader.h:65
itk::ImageSeriesReader::DictionaryArrayType
std::vector< DictionaryRawPointer > DictionaryArrayType
Definition: itkImageSeriesReader.h:76
itk::ImageSeriesReader::FileNamesContainer
std::vector< std::string > FileNamesContainer
Definition: itkImageSeriesReader.h:79
itk::TimeStamp
Generate a unique, increasing time value.
Definition: itkTimeStamp.h:60
itk::ImageSeriesReader::DictionaryArrayRawPointer
const DictionaryArrayType * DictionaryArrayRawPointer
Definition: itkImageSeriesReader.h:77
itk::ImageSeriesReader::AddFileName
void AddFileName(std::string const &name)
Definition: itkImageSeriesReader.h:113
itkMetaDataDictionary.h
itk::ImageSeriesReader::m_FileNames
FileNamesContainer m_FileNames
Definition: itkImageSeriesReader.h:207
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itk::ImageSource::OutputImagePixelType
typename OutputImageType::PixelType OutputImagePixelType
Definition: itkImageSource.h:93
itk::Math::e
static constexpr double e
Definition: itkMath.h:54
itkSize.h
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293