ITK  5.2.0
Insight Toolkit
itkNiftiImageIO.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 
19 #ifndef itkNiftiImageIO_h
20 #define itkNiftiImageIO_h
21 #include "ITKIONIFTIExport.h"
22 
23 
24 #include <fstream>
25 #include <memory>
26 #include "itkImageIOBase.h"
27 
28 namespace itk
29 {
30 
36 {
37 public:
42  enum class Analyze75Flavor : uint8_t
43  {
45  AnalyzeITK4 = 4,
46 
48  AnalyzeFSL = 3,
49 
51  AnalyzeSPM = 2,
52 
55 
57  AnalyzeReject = 0
58  };
59 
64  enum class NiftiFileEnum : int8_t
65  {
67  TwoFileNifti = 2,
68 
70  OneFileNifti = 1,
71 
73  Analyze75 = 0,
74 
76  OtherOrError = -1,
77  };
78 };
79 
81 #if !defined(ITK_LEGACY_REMOVE)
83 #endif
84 
86 extern ITKIONIFTI_EXPORT std::ostream &
87  operator<<(std::ostream & out, const NiftiImageIOEnums::Analyze75Flavor value);
88 extern ITKIONIFTI_EXPORT std::ostream &
89  operator<<(std::ostream & out, const NiftiImageIOEnums::NiftiFileEnum value);
91 
105 class ITKIONIFTI_EXPORT NiftiImageIO : public ImageIOBase
106 {
107 public:
108  ITK_DISALLOW_COPY_AND_MOVE(NiftiImageIO);
109 
114 
116  itkNewMacro(Self);
117 
119  itkTypeMacro(NiftiImageIO, Superclass);
120 
121  //-------- This part of the interfaces deals with reading data. -----
122 
123 #if !defined(ITK_LEGACY_REMOVE)
124 
125  using FileType = NiftiImageIOEnums::NiftiFileEnum;
126  // We need to expose the enum values at the class level
127  // for backwards compatibility
128  static constexpr FileType TwoFileNifti = NiftiImageIOEnums::NiftiFileEnum::TwoFileNifti;
129  static constexpr FileType OneFileNifti = NiftiImageIOEnums::NiftiFileEnum::OneFileNifti;
130  static constexpr FileType Analyze75 = NiftiImageIOEnums::NiftiFileEnum::Analyze75;
131  static constexpr FileType OtherOrError = NiftiImageIOEnums::NiftiFileEnum::OtherOrError;
132 #endif
133 
141  DetermineFileType(const char * FileNameToRead);
142 
149  bool
150  CanReadFile(const char * FileNameToRead) override;
151 
153  void
154  ReadImageInformation() override;
155 
157  void
158  Read(void * buffer) override;
159 
160  //-------- This part of the interfaces deals with writing data. -----
161 
167  bool
168  CanWriteFile(const char * FileNameToWrite) override;
169 
175  void
176  WriteImageInformation() override;
177 
180  void
181  Write(const void * buffer) override;
182 
186  GenerateStreamableReadRegionFromRequestedRegion(const ImageIORegion & requestedRegion) const override;
187 
189  itkSetMacro(RescaleSlope, double);
190  itkSetMacro(RescaleIntercept, double);
192 
197  itkSetMacro(LegacyAnalyze75Mode, NiftiImageIOEnums::Analyze75Flavor);
198  itkGetConstMacro(LegacyAnalyze75Mode, NiftiImageIOEnums::Analyze75Flavor);
200 
201 protected:
202  NiftiImageIO();
203  ~NiftiImageIO() override;
204  void
205  PrintSelf(std::ostream & os, Indent indent) const override;
206 
207  virtual bool
209  {
210  return false;
211  }
212 
213 private:
214  // Try to use the Q and S form codes from MetaDataDictionary if they are specified
215  // there, otherwise default to the backwards compatible values from earlier
216  // versions of ITK. The qform guess would probably been better to have
217  // been guessed as NIFTI_XFORM_SCANNER_ANAT
218  unsigned int
219  getSFormCodeFromDictionary() const;
220  unsigned int
221  getQFormCodeFromDictionary() const;
222 
223  bool
224  MustRescale() const;
225 
226  void
227  DefineHeaderObjectDataType();
228 
229  void
230  SetNIfTIOrientationFromImageIO(unsigned short int origdims, unsigned short int dims);
231 
232  void
233  SetImageIOOrientationFromNIfTI(unsigned short int dims);
234 
235  void
236  SetImageIOMetadataFromNIfTI();
237 
238  // This proxy class provides a nifti_image pointer interface to the internal implementation
239  // of itk::NiftiImageIO, while hiding the niftilib interface from the external ITK interface.
240  class NiftiImageProxy;
241 
242  // Note that it is essential that m_NiftiImageHolder is defined before m_NiftiImage, to ensure that
243  // m_NiftiImage can directly get a proxy from m_NiftiImageHolder during NiftiImageIO construction.
244  const std::unique_ptr<NiftiImageProxy> m_NiftiImageHolder;
245 
246  NiftiImageProxy & m_NiftiImage;
247 
248  double m_RescaleSlope{ 1.0 };
249  double m_RescaleIntercept{ 0.0 };
250 
252 
254 };
255 
256 
257 } // end namespace itk
258 
259 #endif // itkNiftiImageIO_h
itk::ImageIOBase
Abstract superclass defines image IO interface.
Definition: itkImageIOBase.h:77
itk::uint8_t
::uint8_t uint8_t
Definition: itkIntTypes.h:29
itk::NiftiImageIOEnums::NiftiFileEnum::OneFileNifti
itk::NiftiImageIOEnums::NiftiFileEnum::TwoFileNifti
itk::NiftiImageIOEnums::Analyze75Flavor::AnalyzeFSL
itk::operator<<
std::ostream & operator<<(std::ostream &os, const Array< TValue > &arr)
Definition: itkArray.h:218
itk::int8_t
::int8_t int8_t
Definition: itkIntTypes.h:28
itk::NiftiImageIO::m_NiftiImage
NiftiImageProxy & m_NiftiImage
Definition: itkNiftiImageIO.h:246
itk::NiftiImageIOEnums::NiftiFileEnum
NiftiFileEnum
Definition: itkNiftiImageIO.h:64
itk::SmartPointer< Self >
itk::CommonEnums::IOComponent
IOComponent
Definition: itkCommonEnums.h:76
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::NiftiImageIO::m_NiftiImageHolder
const std::unique_ptr< NiftiImageProxy > m_NiftiImageHolder
Definition: itkNiftiImageIO.h:240
itk::ImageIORegion
An ImageIORegion represents a structured region of data.
Definition: itkImageIORegion.h:52
itk::NiftiImageIOEnums::Analyze75Flavor::AnalyzeReject
itk::NiftiImageIO::GetUseLegacyModeForTwoFileWriting
virtual bool GetUseLegacyModeForTwoFileWriting() const
Definition: itkNiftiImageIO.h:208
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itk::NiftiImageIOEnums::NiftiFileEnum::OtherOrError
itk::NiftiImageIOEnums::NiftiFileEnum::Analyze75
itkImageIOBase.h
NiftiFileEnum
itk::NiftiImageIOEnums::Analyze75Flavor
Analyze75Flavor
Definition: itkNiftiImageIO.h:42
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::NiftiImageIOEnums::Analyze75Flavor::AnalyzeSPM
itk::NiftiImageIOEnums::Analyze75Flavor::AnalyzeITK4Warning
itk::NiftiImageIO
Class that defines how to read Nifti file format. Nifti IMAGE FILE FORMAT - As much information as I ...
Definition: itkNiftiImageIO.h:105
itk::NiftiImageIOEnums::Analyze75Flavor::AnalyzeITK4
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:62
itk::CommonEnums::IOComponent::UNKNOWNCOMPONENTTYPE
Analyze75Flavor
itk::NiftiImageIOEnums
Definition: itkNiftiImageIO.h:35
itk::NiftiImageIO::m_LegacyAnalyze75Mode
NiftiImageIOEnums::Analyze75Flavor m_LegacyAnalyze75Mode
Definition: itkNiftiImageIO.h:253