ITK  6.0.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 * 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
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
28namespace itk
29{
30
36{
37public:
42 enum class Analyze75Flavor : uint8_t
43 {
44
46 AnalyzeITK4 = 4,
47
49 AnalyzeFSL = 3,
50
52 AnalyzeSPM = 2,
53
55 AnalyzeITK4Warning = 1,
56
58 AnalyzeReject = 0
59 };
60
65 enum class NiftiFileEnum : int8_t
66 {
68 TwoFileNifti = 2,
69
71 OneFileNifti = 1,
72
74 Analyze75 = 0,
75
77 OtherOrError = -1,
78 };
79};
80
82#if !defined(ITK_LEGACY_REMOVE)
83using Analyze75Flavor = NiftiImageIOEnums::Analyze75Flavor;
84#endif
85
87extern ITKIONIFTI_EXPORT std::ostream &
88 operator<<(std::ostream & out, const NiftiImageIOEnums::Analyze75Flavor value);
89extern ITKIONIFTI_EXPORT std::ostream &
90 operator<<(std::ostream & out, const NiftiImageIOEnums::NiftiFileEnum value);
106class ITKIONIFTI_EXPORT NiftiImageIO : public ImageIOBase
107{
108public:
109 ITK_DISALLOW_COPY_AND_MOVE(NiftiImageIO);
110
115
117 itkNewMacro(Self);
118
120 itkOverrideGetNameOfClassMacro(NiftiImageIO);
121
122 //-------- This part of the interfaces deals with reading data. -----
123
124#if !defined(ITK_LEGACY_REMOVE)
126 using FileType = NiftiImageIOEnums::NiftiFileEnum;
127 // We need to expose the enum values at the class level
128 // for backwards compatibility
129 static constexpr FileType TwoFileNifti = NiftiImageIOEnums::NiftiFileEnum::TwoFileNifti;
130 static constexpr FileType OneFileNifti = NiftiImageIOEnums::NiftiFileEnum::OneFileNifti;
131 static constexpr FileType Analyze75 = NiftiImageIOEnums::NiftiFileEnum::Analyze75;
132 static constexpr FileType OtherOrError = NiftiImageIOEnums::NiftiFileEnum::OtherOrError;
133#endif
134
142 DetermineFileType(const char * FileNameToRead);
143
150 bool
151 CanReadFile(const char * FileNameToRead) override;
152
154 void
156
158 void
159 Read(void * buffer) override;
160
161 //-------- This part of the interfaces deals with writing data. -----
162
168 bool
169 CanWriteFile(const char * FileNameToWrite) override;
170
176 void
178
181 void
182 Write(const void * buffer) override;
183
187 GenerateStreamableReadRegionFromRequestedRegion(const ImageIORegion & requestedRegion) const override;
188
190 itkSetMacro(RescaleSlope, double);
191 itkSetMacro(RescaleIntercept, double);
198 itkSetMacro(LegacyAnalyze75Mode, NiftiImageIOEnums::Analyze75Flavor);
199 itkGetConstMacro(LegacyAnalyze75Mode, NiftiImageIOEnums::Analyze75Flavor);
210 itkSetMacro(ConvertRASVectors, bool);
211 itkGetConstMacro(ConvertRASVectors, bool);
212 itkBooleanMacro(ConvertRASVectors);
230 itkSetMacro(ConvertRASDisplacementVectors, bool);
231 itkGetConstMacro(ConvertRASDisplacementVectors, bool);
232 itkBooleanMacro(ConvertRASDisplacementVectors);
236 itkSetMacro(SFORM_Permissive, bool);
237 itkGetConstMacro(SFORM_Permissive, bool);
238 itkBooleanMacro(SFORM_Permissive);
241protected:
243 ~NiftiImageIO() override;
244 void
245 PrintSelf(std::ostream & os, Indent indent) const override;
246
247 virtual bool
249 {
250 return false;
251 }
252
253private:
254 // Try to use the Q and S form codes from MetaDataDictionary if they are specified
255 // there, otherwise default to the backwards compatible values from earlier
256 // versions of ITK. The qform guess would probably been better to have
257 // been guessed as NIFTI_XFORM_SCANNER_ANAT
258 unsigned int
260 unsigned int
262
263 bool
264 MustRescale() const;
265
266 void
268
269 void
270 SetNIfTIOrientationFromImageIO(unsigned short origdims, unsigned short dims);
271
272 void
273 SetImageIOOrientationFromNIfTI(unsigned short dims, double spacingscale, double timingscale);
274
275 void
277
278 // This proxy class provides a nifti_image pointer interface to the internal implementation
279 // of itk::NiftiImageIO, while hiding the niftilib interface from the external ITK interface.
280 class NiftiImageProxy;
281
282 // Note that it is essential that m_NiftiImageHolder is defined before m_NiftiImage, to ensure that
283 // m_NiftiImage can directly get a proxy from m_NiftiImageHolder during NiftiImageIO construction.
284 const std::unique_ptr<NiftiImageProxy> m_NiftiImageHolder;
285
286 NiftiImageProxy & m_NiftiImage;
287
288 double m_RescaleSlope{ 1.0 };
289 double m_RescaleIntercept{ 0.0 };
290
291 bool m_ConvertRAS{ false };
292 bool m_ConvertRASVectors{ false };
293 bool m_ConvertRASDisplacementVectors{ true };
294
296
297 NiftiImageIOEnums::Analyze75Flavor m_LegacyAnalyze75Mode{};
298
300 bool m_SFORM_Corrected{ false };
301};
302
303
304} // end namespace itk
305
306#endif // itkNiftiImageIO_h
Abstract superclass defines image IO interface.
An ImageIORegion represents a structured region of data.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
Class that defines how to read Nifti file format. Nifti IMAGE FILE FORMAT - As much information as I ...
ImageIORegion GenerateStreamableReadRegionFromRequestedRegion(const ImageIORegion &requestedRegion) const override
void SetImageIOOrientationFromNIfTI(unsigned short dims, double spacingscale, double timingscale)
void ReadImageInformation() override
NiftiImageProxy & m_NiftiImage
void SetImageIOMetadataFromNIfTI()
const std::unique_ptr< NiftiImageProxy > m_NiftiImageHolder
void PrintSelf(std::ostream &os, Indent indent) const override
void SetNIfTIOrientationFromImageIO(unsigned short origdims, unsigned short dims)
unsigned int getQFormCodeFromDictionary() const
virtual bool GetUseLegacyModeForTwoFileWriting() const
bool CanWriteFile(const char *FileNameToWrite) override
void WriteImageInformation() override
NiftiImageIOEnums::NiftiFileEnum DetermineFileType(const char *FileNameToRead)
void Read(void *buffer) override
bool MustRescale() const
~NiftiImageIO() override
void Write(const void *buffer) override
unsigned int getSFormCodeFromDictionary() const
void DefineHeaderObjectDataType()
bool CanReadFile(const char *FileNameToRead) override
Base class for most ITK classes.
Definition: itkObject.h:62
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)