ITK  5.4.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 {
45 AnalyzeITK4 = 4,
46
48 AnalyzeFSL = 3,
49
51 AnalyzeSPM = 2,
52
54 AnalyzeITK4Warning = 1,
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
86extern ITKIONIFTI_EXPORT std::ostream &
87 operator<<(std::ostream & out, const NiftiImageIOEnums::Analyze75Flavor value);
88extern ITKIONIFTI_EXPORT std::ostream &
89 operator<<(std::ostream & out, const NiftiImageIOEnums::NiftiFileEnum value);
105class ITKIONIFTI_EXPORT NiftiImageIO : public ImageIOBase
106{
107public:
108 ITK_DISALLOW_COPY_AND_MOVE(NiftiImageIO);
109
114
116 itkNewMacro(Self);
117
119 itkOverrideGetNameOfClassMacro(NiftiImageIO);
120
121 //-------- This part of the interfaces deals with reading data. -----
122
123#if !defined(ITK_LEGACY_REMOVE)
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
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
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);
197 itkSetMacro(LegacyAnalyze75Mode, NiftiImageIOEnums::Analyze75Flavor);
198 itkGetConstMacro(LegacyAnalyze75Mode, NiftiImageIOEnums::Analyze75Flavor);
209 itkSetMacro(ConvertRASVectors, bool);
210 itkGetConstMacro(ConvertRASVectors, bool);
211 itkBooleanMacro(ConvertRASVectors);
229 itkSetMacro(ConvertRASDisplacementVectors, bool);
230 itkGetConstMacro(ConvertRASDisplacementVectors, bool);
231 itkBooleanMacro(ConvertRASDisplacementVectors);
235 itkSetMacro(SFORM_Permissive, bool);
236 itkGetConstMacro(SFORM_Permissive, bool);
237 itkBooleanMacro(SFORM_Permissive);
240protected:
242 ~NiftiImageIO() override;
243 void
244 PrintSelf(std::ostream & os, Indent indent) const override;
245
246 virtual bool
248 {
249 return false;
250 }
251
252private:
253 // Try to use the Q and S form codes from MetaDataDictionary if they are specified
254 // there, otherwise default to the backwards compatible values from earlier
255 // versions of ITK. The qform guess would probably been better to have
256 // been guessed as NIFTI_XFORM_SCANNER_ANAT
257 unsigned int
259 unsigned int
261
262 bool
263 MustRescale() const;
264
265 void
267
268 void
269 SetNIfTIOrientationFromImageIO(unsigned short origdims, unsigned short dims);
270
271 void
272 SetImageIOOrientationFromNIfTI(unsigned short dims, double spacingscale, double timingscale);
273
274 void
276
277 // This proxy class provides a nifti_image pointer interface to the internal implementation
278 // of itk::NiftiImageIO, while hiding the niftilib interface from the external ITK interface.
279 class NiftiImageProxy;
280
281 // Note that it is essential that m_NiftiImageHolder is defined before m_NiftiImage, to ensure that
282 // m_NiftiImage can directly get a proxy from m_NiftiImageHolder during NiftiImageIO construction.
283 const std::unique_ptr<NiftiImageProxy> m_NiftiImageHolder;
284
285 NiftiImageProxy & m_NiftiImage;
286
287 double m_RescaleSlope{ 1.0 };
288 double m_RescaleIntercept{ 0.0 };
289
290 bool m_ConvertRAS{ false };
291 bool m_ConvertRASVectors{ false };
292 bool m_ConvertRASDisplacementVectors{ true };
293
295
296 NiftiImageIOEnums::Analyze75Flavor m_LegacyAnalyze75Mode{};
297
299 bool m_SFORM_Corrected{ false };
300};
301
302
303} // end namespace itk
304
305#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....
std::ostream & operator<<(std::ostream &os, const Array< TValue > &arr)
Definition: itkArray.h:216