ITK  5.2.0
Insight Toolkit
itkRawImageIO.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 itkRawImageIO_h
19 #define itkRawImageIO_h
20 
21 #include "itkImageIOBase.h"
22 #include "itkImageRegion.h"
23 #include "itkPixelTraits.h"
24 #include "itkByteSwapper.h"
25 #include "itkVersion.h"
26 #include <string>
27 #include <fstream>
28 
29 namespace itk
30 {
48 template <typename TPixel, unsigned int VImageDimension = 2>
49 class ITK_TEMPLATE_EXPORT RawImageIO : public ImageIOBase
50 {
51 public:
52  ITK_DISALLOW_COPY_AND_MOVE(RawImageIO);
53 
55  using Self = RawImageIO;
58 
60  itkNewMacro(Self);
61 
63  itkTypeMacro(RawImageIO, ImageIOBase);
64 
67  using PixelType = TPixel;
68 
71 
74 
77 
80  void
81  SetHeaderSize(SizeValueType size);
82 
84  GetHeaderSize();
85 
89  itkSetMacro(FileDimensionality, unsigned long);
90  itkGetConstMacro(FileDimensionality, unsigned long);
92 
98  bool
99  SupportsDimension(unsigned long dim) override
100  {
101  return (dim == m_FileDimensionality);
102  }
103 
104  /*-------- This part of the interface deals with reading data. ------ */
105 
109  bool
110  CanReadFile(const char *) override
111  {
112  return false;
113  }
114 
117  void
119  {
120  return;
121  }
122 
124  void
125  Read(void * buffer) override;
126 
128  itkGetConstReferenceMacro(ImageMask, unsigned short);
129  void
130  SetImageMask(unsigned long val)
131  {
132  if (val == m_ImageMask)
133  {
134  return;
135  }
136  m_ImageMask = ((unsigned short)(val));
137  this->Modified();
138  }
140 
142  virtual void
143  ReadHeader(const std::string = std::string())
144  {}
145 
146  /*-------- This part of the interfaces deals with writing data. ----- */
147 
151  bool
152  CanWriteFile(const char *) override;
153 
155  void
157  {
158  return;
159  }
160 
162  void
163  Write(const void * buffer) override;
164 
165 protected:
166  RawImageIO();
167  ~RawImageIO() override = default;
168  void
169  PrintSelf(std::ostream & os, Indent indent) const override;
170 
171  // void ComputeInternalFileName(unsigned long slice);
172 
173 private:
174  std::string m_InternalFileName;
175 
176  unsigned long m_FileDimensionality;
179  unsigned short m_ImageMask;
180 };
181 
182 template <typename TPixel, unsigned int VImageDimension>
183 class ITK_TEMPLATE_EXPORT RawImageIOFactory : public ObjectFactoryBase
184 {
185 public:
186  ITK_DISALLOW_COPY_AND_MOVE(RawImageIOFactory);
187 
193 
195  const char *
196  GetITKSourceVersion() const override
197  {
198  return ITK_SOURCE_VERSION;
199  }
200 
201  const char *
202  GetDescription() const override
203  {
204  return "Raw ImageIO Factory, allows the loading of Raw images into insight";
205  }
206 
208  itkFactorylessNewMacro(Self);
209 
211  itkTypeMacro(RawImageIOFactory, ObjectFactoryBase);
212 
214  static void
216  {
218  }
219 
220 protected:
221  RawImageIOFactory() = default;
222  ~RawImageIOFactory() override = default;
223 };
224 } // namespace itk
225 
226 #ifndef ITK_MANUAL_INSTANTIATION
227 # include "itkRawImageIO.hxx"
228 #endif
229 
230 #endif
itk::ImageIOBase
Abstract superclass defines image IO interface.
Definition: itkImageIOBase.h:77
itk::ObjectFactoryBase
Create instances of classes using an object factory.
Definition: itkObjectFactoryBase.h:62
itkByteSwapper.h
itk::RawImageIO::WriteImageInformation
void WriteImageInformation() override
Definition: itkRawImageIO.h:156
itk::RawImageIO::ReadImageInformation
void ReadImageInformation() override
Definition: itkRawImageIO.h:118
itk::RawImageIOFactory
Definition: itkRawImageIO.h:183
itk::RawImageIO::ComponentType
typename PixelTraits< PixelType >::ValueType ComponentType
Definition: itkRawImageIO.h:73
itkPixelTraits.h
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::RawImageIO::m_ImageMask
unsigned short m_ImageMask
Definition: itkRawImageIO.h:179
itk::RawImageIO::PixelType
TPixel PixelType
Definition: itkRawImageIO.h:67
itk::RawImageIO::m_FileDimensionality
unsigned long m_FileDimensionality
Definition: itkRawImageIO.h:176
itkImageRegion.h
ITK_SOURCE_VERSION
#define ITK_SOURCE_VERSION
Definition: itkVersion.h:39
itk::RawImageIO::CanReadFile
bool CanReadFile(const char *) override
Definition: itkRawImageIO.h:110
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itk::RawImageIOFactory::GetITKSourceVersion
const char * GetITKSourceVersion() const override
Definition: itkRawImageIO.h:196
itk::RawImageIOFactory::GetDescription
const char * GetDescription() const override
Definition: itkRawImageIO.h:202
itk::ImageIOBase::SizeValueType
::itk::SizeValueType SizeValueType
Definition: itkImageIOBase.h:98
itk::ByteSwapper
Perform machine dependent byte swapping.
Definition: itkByteSwapper.h:50
itk::RawImageIO::ReadHeader
virtual void ReadHeader(const std::string=std::string())
Definition: itkRawImageIO.h:143
itk::RawImageIO::m_ManualHeaderSize
bool m_ManualHeaderSize
Definition: itkRawImageIO.h:177
itkImageIOBase.h
itk::RawImageIO::SetImageMask
void SetImageMask(unsigned long val)
Definition: itkRawImageIO.h:130
itk::RawImageIO::m_InternalFileName
std::string m_InternalFileName
Definition: itkRawImageIO.h:174
itkVersion.h
itk::PixelTraits::ValueType
typename TPixelType::ValueType ValueType
Definition: itkPixelTraits.h:52
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:62
itk::RawImageIO::m_HeaderSize
SizeValueType m_HeaderSize
Definition: itkRawImageIO.h:178
itk::RawImageIO
Read and write raw binary images.
Definition: itkRawImageIO.h:49
itk::ObjectFactoryBase::RegisterFactory
static bool RegisterFactory(ObjectFactoryBase *, InsertionPositionEnum where=InsertionPositionEnum::INSERT_AT_BACK, vcl_size_t position=0)
itk::RawImageIO::SupportsDimension
bool SupportsDimension(unsigned long dim) override
Definition: itkRawImageIO.h:99
itk::RawImageIOFactory::RegisterOneFactory
static void RegisterOneFactory()
Definition: itkRawImageIO.h:215
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83