ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
28
29namespace itk
30{
31
37{
38public:
43 enum class Analyze75Flavor : uint8_t
44 {
45
48
51
54
57
60 };
61
66 enum class NiftiFileEnum : int8_t
67 {
70
73
76
79 };
80};
81
83#if !defined(ITK_LEGACY_REMOVE)
84using Analyze75Flavor = NiftiImageIOEnums::Analyze75Flavor;
85#endif
86
89extern ITKIONIFTI_EXPORT std::ostream &
90 operator<<(std::ostream & out, const NiftiImageIOEnums::Analyze75Flavor value);
91extern ITKIONIFTI_EXPORT std::ostream &
92 operator<<(std::ostream & out, const NiftiImageIOEnums::NiftiFileEnum value);
107class ITKIONIFTI_EXPORT NiftiImageIO : public ImageIOBase
108{
109public:
110 ITK_DISALLOW_COPY_AND_MOVE(NiftiImageIO);
111
116
118 itkNewMacro(Self);
119
121 itkOverrideGetNameOfClassMacro(NiftiImageIO);
122
123 //-------- This part of the interfaces deals with reading data. -----
124
125#if !defined(ITK_LEGACY_REMOVE)
127 using FileType = NiftiImageIOEnums::NiftiFileEnum;
128 // We need to expose the enum values at the class level
129 // for backwards compatibility
130 static constexpr FileType TwoFileNifti = NiftiImageIOEnums::NiftiFileEnum::TwoFileNifti;
131 static constexpr FileType OneFileNifti = NiftiImageIOEnums::NiftiFileEnum::OneFileNifti;
132 static constexpr FileType Analyze75 = NiftiImageIOEnums::NiftiFileEnum::Analyze75;
133 static constexpr FileType OtherOrError = NiftiImageIOEnums::NiftiFileEnum::OtherOrError;
134#endif
135
143 DetermineFileType(const char * FileNameToRead);
144
151 bool
152 CanReadFile(const char * FileNameToRead) override;
153
155 void
157
159 void
160 Read(void * buffer) override;
161
162 //-------- This part of the interfaces deals with writing data. -----
163
169 bool
170 CanWriteFile(const char * FileNameToWrite) override;
171
177 void
179
182 void
183 Write(const void * buffer) override;
184
188 GenerateStreamableReadRegionFromRequestedRegion(const ImageIORegion & requestedRegion) const override;
189
192 itkSetMacro(RescaleSlope, double);
193 itkSetMacro(RescaleIntercept, double);
200 itkSetMacro(LegacyAnalyze75Mode, NiftiImageIOEnums::Analyze75Flavor);
201 itkGetConstMacro(LegacyAnalyze75Mode, NiftiImageIOEnums::Analyze75Flavor);
212 itkSetMacro(ConvertRASVectors, bool);
213 itkGetConstMacro(ConvertRASVectors, bool);
214 itkBooleanMacro(ConvertRASVectors);
232 itkSetMacro(ConvertRASDisplacementVectors, bool);
233 itkGetConstMacro(ConvertRASDisplacementVectors, bool);
234 itkBooleanMacro(ConvertRASDisplacementVectors);
238 itkSetMacro(SFORM_Permissive, bool);
239 itkGetConstMacro(SFORM_Permissive, bool);
240 itkBooleanMacro(SFORM_Permissive);
242protected:
244 ~NiftiImageIO() override;
245 void
246 PrintSelf(std::ostream & os, Indent indent) const override;
247
248 virtual bool
250 {
251 return false;
252 }
253
254private:
255 // This proxy class has a nifti_image as a smart pointer to encapsulate access to the header.
256 class NiftiImageProxy;
257
258 std::unique_ptr<NiftiImageProxy> m_Holder;
259
260 // Try to use the Q and S form codes from MetaDataDictionary if they are specified
261 // there, otherwise default to the backwards compatible values from earlier
262 // versions of ITK. The qform guess would probably been better to have
263 // been guessed as NIFTI_XFORM_SCANNER_ANAT
264 unsigned int
266 unsigned int
268
269 bool
270 MustRescale() const;
271
272 void
274
275 void
276 SetNIfTIOrientationFromImageIO(unsigned short origdims, unsigned short dims);
277
278 void
279 SetImageIOOrientationFromNIfTI(unsigned short dims, double spacingscale, double timingscale);
280
281 void
283
284 double m_RescaleSlope{ 1.0 };
285 double m_RescaleIntercept{ 0.0 };
286
287 bool m_ConvertRAS{ false };
288 bool m_ConvertRASVectors{ false };
290
291 IOComponentEnum m_OnDiskComponentType{ IOComponentEnum::UNKNOWNCOMPONENTTYPE };
292
294
295 bool m_SFORM_Permissive{ false };
296 bool m_SFORM_Corrected{ false };
297};
298
299
300} // end namespace itk
301
302#endif // itkNiftiImageIO_h
itk::IOComponentEnum IOComponentEnum
An ImageIORegion represents a structured region of data.
Control indentation during Print() invocation.
Definition itkIndent.h:50
ImageIORegion GenerateStreamableReadRegionFromRequestedRegion(const ImageIORegion &requestedRegion) const override
void SetImageIOOrientationFromNIfTI(unsigned short dims, double spacingscale, double timingscale)
IOComponentEnum m_OnDiskComponentType
std::unique_ptr< NiftiImageProxy > m_Holder
void ReadImageInformation() override
void SetImageIOMetadataFromNIfTI()
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
SmartPointer< Self > Pointer
void WriteImageInformation() override
NiftiImageIOEnums::NiftiFileEnum DetermineFileType(const char *FileNameToRead)
void Read(void *buffer) override
bool MustRescale() const
NiftiImageIOEnums::Analyze75Flavor m_LegacyAnalyze75Mode
~NiftiImageIO() override
void Write(const void *buffer) override
unsigned int getSFormCodeFromDictionary() const
void DefineHeaderObjectDataType()
bool CanReadFile(const char *FileNameToRead) override
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)