ITK  5.1.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 
31 
37 {
39  AnalyzeITK4 = 4,
40 
42  AnalyzeFSL = 3,
43 
45  AnalyzeSPM = 2,
46 
49 
51  AnalyzeReject = 0
52 };
53 
55 extern ITKIONIFTI_EXPORT std::ostream &
56  operator<<(std::ostream & out, const Analyze75Flavor value);
57 
71 class ITKIONIFTI_EXPORT NiftiImageIO : public ImageIOBase
72 {
73 public:
74  ITK_DISALLOW_COPY_AND_ASSIGN(NiftiImageIO);
75 
77  using Self = NiftiImageIO;
80 
82  itkNewMacro(Self);
83 
85  itkTypeMacro(NiftiImageIO, Superclass);
86 
87  //-------- This part of the interfaces deals with reading data. -----
88 
91  enum FileType
92  {
93 
95  TwoFileNifti = 2,
96 
98  OneFileNifti = 1,
99 
101  Analyze75 = 0,
102 
104  OtherOrError = -1,
105  };
106 
107 
114  FileType
115  DetermineFileType(const char * FileNameToRead);
116 
123  bool
124  CanReadFile(const char * FileNameToRead) override;
125 
127  void
128  ReadImageInformation() override;
129 
131  void
132  Read(void * buffer) override;
133 
134  //-------- This part of the interfaces deals with writing data. -----
135 
141  bool
142  CanWriteFile(const char * FileNameToWrite) override;
143 
149  void
150  WriteImageInformation() override;
151 
154  void
155  Write(const void * buffer) override;
156 
160  GenerateStreamableReadRegionFromRequestedRegion(const ImageIORegion & requestedRegion) const override;
161 
163  itkSetMacro(RescaleSlope, double);
164  itkSetMacro(RescaleIntercept, double);
166 
171  itkSetMacro(LegacyAnalyze75Mode, Analyze75Flavor);
172  itkGetConstMacro(LegacyAnalyze75Mode, Analyze75Flavor);
174 
175 protected:
176  NiftiImageIO();
177  ~NiftiImageIO() override;
178  void
179  PrintSelf(std::ostream & os, Indent indent) const override;
180 
181  virtual bool
183  {
184  return false;
185  }
186 
187 private:
188  // Try to use the Q and S form codes from MetaDataDictionary if they are specified
189  // there, otherwise default to the backwards compatible values from earlier
190  // versions of ITK. The qform guess would probably been better to have
191  // been guessed as NIFTI_XFORM_SCANNER_ANAT
192  unsigned int
193  getSFormCodeFromDictionary() const;
194  unsigned int
195  getQFormCodeFromDictionary() const;
196 
197  bool
198  MustRescale() const;
199 
200  void
201  DefineHeaderObjectDataType();
202 
203  void
204  SetNIfTIOrientationFromImageIO(unsigned short int origdims, unsigned short int dims);
205 
206  void
207  SetImageIOOrientationFromNIfTI(unsigned short int dims);
208 
209  void
210  SetImageIOMetadataFromNIfTI();
211 
212  // This proxy class provides a nifti_image pointer interface to the internal implementation
213  // of itk::NiftiImageIO, while hiding the niftilib interface from the external ITK interface.
214  class NiftiImageProxy;
215 
216  // Note that it is essential that m_NiftiImageHolder is defined before m_NiftiImage, to ensure that
217  // m_NiftiImage can directly get a proxy from m_NiftiImageHolder during NiftiImageIO construction.
218  const std::unique_ptr<NiftiImageProxy> m_NiftiImageHolder;
219 
220  NiftiImageProxy & m_NiftiImage;
221 
222  double m_RescaleSlope{ 1.0 };
223  double m_RescaleIntercept{ 0.0 };
224 
226 
228 };
229 
230 
231 } // end namespace itk
232 
233 #endif // itkNiftiImageIO_h
itk::ImageIOBase
Abstract superclass defines image IO interface.
Definition: itkImageIOBase.h:75
itk::uint8_t
::uint8_t uint8_t
Definition: itkIntTypes.h:29
itk::operator<<
std::ostream & operator<<(std::ostream &os, const Array< TValue > &arr)
Definition: itkArray.h:213
itk::Analyze75Flavor::AnalyzeReject
itk::NiftiImageIO::m_NiftiImage
NiftiImageProxy & m_NiftiImage
Definition: itkNiftiImageIO.h:220
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:214
itk::ImageIORegion
An ImageIORegion represents a structured region of data.
Definition: itkImageIORegion.h:52
itk::NiftiImageIO::GetUseLegacyModeForTwoFileWriting
virtual bool GetUseLegacyModeForTwoFileWriting() const
Definition: itkNiftiImageIO.h:182
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itk::NiftiImageIO::m_LegacyAnalyze75Mode
Analyze75Flavor m_LegacyAnalyze75Mode
Definition: itkNiftiImageIO.h:227
itk::Analyze75Flavor::AnalyzeSPM
itkImageIOBase.h
itk::Analyze75Flavor::AnalyzeITK4
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkArray.h:26
itk::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:71
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:60
itk::CommonEnums::IOComponent::UNKNOWNCOMPONENTTYPE
Analyze75Flavor
itk::NiftiImageIO::FileType
FileType
Definition: itkNiftiImageIO.h:91
itk::Analyze75Flavor::AnalyzeFSL