ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkImageFileReader.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 itkImageFileReader_h
19#define itkImageFileReader_h
21
22#include "ITKIOImageBaseExport.h"
23
24#include "itkImageIOBase.h"
25#include "itkImageSource.h"
26#include "itkMacro.h"
27#include "itkImageRegion.h"
30
31namespace itk
32{
33
73template <typename TOutputImage,
75class ITK_TEMPLATE_EXPORT ImageFileReader : public ImageSource<TOutputImage>
76{
77public:
78 ITK_DISALLOW_COPY_AND_MOVE(ImageFileReader);
79
84
86 itkNewMacro(Self);
87
89 itkOverrideGetNameOfClassMacro(ImageFileReader);
90
92 using SizeType = typename TOutputImage::SizeType;
93
95 using IndexType = typename TOutputImage::IndexType;
96
98 using ImageRegionType = typename TOutputImage::RegionType;
99
101 using OutputImagePixelType = typename TOutputImage::InternalPixelType;
102
104 itkSetGetDecoratedInputMacro(FileName, std::string);
105
113 void
115 itkGetModifiableObjectMacro(ImageIO, ImageIOBase);
117
120 itkSetMacro(UseStreaming, bool);
121 itkGetConstReferenceMacro(UseStreaming, bool);
122 itkBooleanMacro(UseStreaming);
124protected:
126 ~ImageFileReader() override = default;
127 void
128 PrintSelf(std::ostream & os, Indent indent) const override;
129
131 void
132 DoConvertBuffer(const void * inputData, size_t numberOfPixels);
133
139 void
141
144 void
146
152 void
154
156 void
157 GenerateData() override;
158
160
161 bool m_UserSpecifiedImageIO{}; // keep track whether the
162 // ImageIO is user specified
163
165
166private:
167 std::string m_ExceptionMessage{};
168
169 // The region that the ImageIO class will return when we ask to
170 // produce the requested region.
172};
173
174
184template <typename TOutputImage,
185 typename ConvertPixelTraits = DefaultConvertPixelTraits<typename TOutputImage::IOPixelType>>
186typename TOutputImage::Pointer
187ReadImage(const std::string & filename)
188{
190 reader->SetFileName(filename);
191 reader->Update();
192 return reader->GetOutput();
193}
194
195} // namespace itk
196
197#ifndef ITK_MANUAL_INSTANTIATION
198# include "itkImageFileReader.hxx"
199#endif
200
201#if defined ITK_IMAGEIO_FACTORY_REGISTER_MANAGER || defined ITK_IO_FACTORY_REGISTER_MANAGER
202# include "itkImageIOFactoryRegisterManager.h"
203#endif
204
205#endif // itkImageFileReader_h
Base class for all data objects in ITK.
Traits class used to by ConvertPixels to convert blocks of pixels.
void EnlargeOutputRequestedRegion(DataObject *output) override
void SetImageIO(ImageIOBase *imageIO)
void GenerateOutputInformation() override
void DoConvertBuffer(const void *inputData, vcl_size_t numberOfPixels)
typename TOutputImage::RegionType ImageRegionType
void GenerateData() override
typename TOutputImage::IndexType IndexType
void TestFileExistanceAndReadability()
typename TOutputImage::InternalPixelType OutputImagePixelType
typename TOutputImage::SizeType SizeType
ImageSource< TOutputImage > Superclass
static Pointer New()
SmartPointer< Self > Pointer
~ImageFileReader() override=default
ImageIOBase::Pointer m_ImageIO
void PrintSelf(std::ostream &os, Indent indent) const override
Abstract superclass defines image IO interface.
SmartPointer< Self > Pointer
An ImageIORegion represents a structured region of data.
Control indentation during Print() invocation.
Definition itkIndent.h:50
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
TOutputImage::Pointer ReadImage(const std::string &filename)