ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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 * https://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>
26#include "itkImageFileReader.h"
27
28namespace itk
29{
43
44template <typename TOutputImage>
45class ITK_TEMPLATE_EXPORT ImageSeriesReader : public ImageSource<TOutputImage>
46{
47public:
48 ITK_DISALLOW_COPY_AND_MOVE(ImageSeriesReader);
49
54
56 itkNewMacro(Self);
57
59 itkOverrideGetNameOfClassMacro(ImageSeriesReader);
60
62 using SizeType = typename TOutputImage::SizeType;
63
65 using IndexType = typename TOutputImage::IndexType;
66
68 using ImageRegionType = typename TOutputImage::RegionType;
69
71 using OutputImagePixelType = typename TOutputImage::PixelType;
72
76 using DictionaryArrayType = std::vector<DictionaryRawPointer>;
78
79 using FileNamesContainer = std::vector<std::string>;
80
84 void
86 {
87 if (m_FileNames != name)
88 {
89 m_FileNames = name;
90 this->Modified();
91 }
92 }
93
94 const FileNamesContainer &
96 {
97 return m_FileNames;
98 }
99
100
103 void
104 SetFileName(const std::string & name)
105 {
106 m_FileNames.clear();
107 m_FileNames.push_back(name);
108 this->Modified();
109 }
110
112 void
113 AddFileName(const std::string & name)
114 {
115 m_FileNames.push_back(name);
116 this->Modified();
117 }
118
122 itkSetMacro(ReverseOrder, bool);
123 itkGetConstMacro(ReverseOrder, bool);
124 itkBooleanMacro(ReverseOrder);
126
130 itkSetMacro(ForceOrthogonalDirection, bool);
131 itkGetConstMacro(ForceOrthogonalDirection, bool);
132 itkBooleanMacro(ForceOrthogonalDirection);
134
140 itkSetObjectMacro(ImageIO, ImageIOBase);
141 itkGetModifiableObjectMacro(ImageIO, ImageIOBase);
143
153 itkSetMacro(MetaDataDictionaryArrayUpdate, bool);
154 itkGetConstMacro(MetaDataDictionaryArrayUpdate, bool);
155 itkBooleanMacro(MetaDataDictionaryArrayUpdate);
156
159 void
161
167 void
169
174
177 itkSetMacro(UseStreaming, bool);
178 itkGetConstReferenceMacro(UseStreaming, bool);
179 itkBooleanMacro(UseStreaming);
181
184 itkSetMacro(SpacingWarningRelThreshold, double);
185 itkGetConstMacro(SpacingWarningRelThreshold, double);
187protected:
189 : m_ImageIO(nullptr)
190
191 {}
192 ~ImageSeriesReader() override = default;
193 void
194 PrintSelf(std::ostream & os, Indent indent) const override;
195
197 void
198 GenerateData() override;
199
202
204 bool m_ReverseOrder{ false };
205
208
211
218
219 bool m_UseStreaming{ true };
220
221 bool m_SpacingDefined{ false };
222
224
225private:
227
228 int
230
234
236 std::vector<MetaDataDictionary> m_InternalMetaDataDictionaries{};
237
240
243};
244} // namespace itk
245
246#ifndef ITK_MANUAL_INSTANTIATION
247# include "itkImageSeriesReader.hxx"
248#endif
249
250#endif // itkImageSeriesReader_h
Base class for all data objects in ITK.
Data source that reads image data from a single file.
Abstract superclass defines image IO interface.
SmartPointer< Self > Pointer
std::vector< std::string > FileNamesContainer
~ImageSeriesReader() override=default
ImageFileReader< TOutputImage > ReaderType
typename TOutputImage::PixelType OutputImagePixelType
std::vector< DictionaryRawPointer > DictionaryArrayType
void GenerateData() override
ImageSource< TOutputImage > Superclass
FileNamesContainer m_FileNames
typename TOutputImage::RegionType ImageRegionType
const DictionaryArrayType * DictionaryArrayRawPointer
void GenerateOutputInformation() override
std::vector< MetaDataDictionary > m_InternalMetaDataDictionaries
typename TOutputImage::SizeType SizeType
MetaDataDictionary * DictionaryRawPointer
void AddFileName(const std::string &name)
typename TOutputImage::IndexType IndexType
MetaDataDictionary DictionaryType
void SetFileName(const std::string &name)
DictionaryArrayRawPointer GetMetaDataDictionaryArray() const
DictionaryArrayType m_MetaDataDictionaryArray
void PrintSelf(std::ostream &os, Indent indent) const override
int ComputeMovingDimensionIndex(ReaderType *reader)
SmartPointer< Self > Pointer
void SetFileNames(const FileNamesContainer &name)
ImageIOBase::Pointer m_ImageIO
void EnlargeOutputRequestedRegion(DataObject *output) override
const FileNamesContainer & GetFileNames() const
Control indentation during Print() invocation.
Definition itkIndent.h:50
Provides a mechanism for storing a collection of arbitrary data types.
virtual void Modified() const
Implements transparent reference counting.
Generate a unique, increasing time value.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....