ITK  6.0.0
Insight Toolkit
itkGDCMImageIO.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 *
20 * Portions of this file are subject to the VTK Toolkit Version 3 copyright.
21 *
22 * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
23 *
24 * For complete copyright, license and disclaimer of warranty information
25 * please refer to the NOTICE file at the top of the ITK source tree.
26 *
27 *=========================================================================*/
28#ifndef itkGDCMImageIO_h
29#define itkGDCMImageIO_h
30
31#include "itkCommonEnums.h"
32#include "itkImageIOBase.h"
33#include "ITKIOGDCMExport.h"
34#include <fstream>
35#include <string>
36
37
38#if !defined(ITK_LEGACY_REMOVE)
39# define ITKIO_DEPRECATED_GDCM1_API
40#endif
41
42
43namespace itk
44{
49{
50public:
55 enum class Compression : uint8_t
56 {
57 JPEG = 0,
58 JPEG2000,
59 JPEGLS,
60 RLE
61 };
62};
63
64// Define how to print enumeration
65extern ITKIOGDCM_EXPORT std::ostream &
66 operator<<(std::ostream & out, const GDCMImageIOEnums::Compression value);
107class InternalHeader;
108class ITKIOGDCM_EXPORT GDCMImageIO : public ImageIOBase
109{
110public:
111 ITK_DISALLOW_COPY_AND_MOVE(GDCMImageIO);
118
120 itkNewMacro(Self);
121
123 itkOverrideGetNameOfClassMacro(GDCMImageIO);
124
125 /*-------- This part of the interface deals with reading data. ------ */
126
129 bool
130 CanReadFile(const char *) override;
131
133 void
135
137 void
138 Read(void * pointer) override;
139
143 itkGetEnumMacro(InternalComponentType, itk::CommonEnums::IOComponent);
144 itkSetEnumMacro(InternalComponentType, itk::CommonEnums::IOComponent);
147 /*-------- This part of the interfaces deals with writing data. ----- */
148
151 bool
152 CanWriteFile(const char *) override;
153
156 void
158
161 void
162 Write(const void * buffer) override;
163
165 itkGetConstMacro(RescaleSlope, double);
166 itkGetConstMacro(RescaleIntercept, double);
172 itkGetStringMacro(UIDPrefix);
173 itkSetStringMacro(UIDPrefix);
177 itkGetStringMacro(StudyInstanceUID);
178 itkGetStringMacro(SeriesInstanceUID);
179 itkGetStringMacro(FrameOfReferenceInstanceUID);
183 itkSetMacro(KeepOriginalUID, bool);
184 itkGetConstMacro(KeepOriginalUID, bool);
185 itkBooleanMacro(KeepOriginalUID);
190 itkSetMacro(LoadPrivateTags, bool);
191 itkGetConstMacro(LoadPrivateTags, bool);
192 itkBooleanMacro(LoadPrivateTags);
197 itkSetMacro(ReadYBRtoRGB, bool);
198 itkGetConstMacro(ReadYBRtoRGB, bool);
199 itkBooleanMacro(ReadYBRtoRGB);
202#if defined(ITKIO_DEPRECATED_GDCM1_API)
208 void
209 GetPatientName(char * name, size_t len = 512);
210
211 void
212 GetPatientID(char * name, size_t len = 512);
213
214 void
215 GetPatientSex(char * name, size_t len = 512);
216
217 void
218 GetPatientAge(char * name, size_t len = 512);
219
220 void
221 GetStudyID(char * name, size_t len = 512);
222
223 void
224 GetPatientDOB(char * name, size_t len = 512);
225
226 void
227 GetStudyDescription(char * name, size_t len = 512);
228
229 void
230 GetBodyPart(char * name, size_t len = 512);
231
232 void
233 GetNumberOfSeriesInStudy(char * name, size_t len = 512);
234
235 void
236 GetNumberOfStudyRelatedSeries(char * name, size_t len = 512);
237
238 void
239 GetStudyDate(char * name, size_t len = 512);
240
241 void
242 GetModality(char * name, size_t len = 512);
243
244 void
245 GetManufacturer(char * name, size_t len = 512);
246
247 void
248 GetInstitution(char * name, size_t len = 512);
249
250 void
251 GetModel(char * name, size_t len = 512);
252
253 void
254 GetScanOptions(char * name, size_t len = 512);
255#endif
256
259 bool
260 GetValueFromTag(const std::string & tag, std::string & value);
261
268 static bool
269 GetLabelFromTag(const std::string & tag, std::string & labelId);
270
271#if defined(ITKIO_DEPRECATED_GDCM1_API)
278 virtual void
279 SetMaxSizeLoadEntry(const long)
280 {}
281
286 virtual void
287 SetLoadSequences(const bool)
288 {}
289 virtual bool
290 GetLoadSequences() const
291 {
292 return true;
293 }
294 virtual void
295 LoadSequencesOn()
296 {}
297 virtual void
298 LoadSequencesOff()
299 {}
309 static void
310 SetLoadSequencesDefault(bool)
311 {}
312 static void
313 LoadSequencesDefaultOn()
314 {}
315 static void
316 LoadSequencesDefaultOff()
317 {}
318 static bool
319 GetLoadSequencesDefault()
320 {
321 return true;
322 }
332 static void
333 SetLoadPrivateTagsDefault(bool)
334 {}
335 static void
336 LoadPrivateTagsDefaultOn()
337 {}
338 static void
339 LoadPrivateTagsDefaultOff()
340 {}
341 static bool
342 GetLoadPrivateTagsDefault()
343 {
344 return true;
345 }
346#endif
351#if !defined(ITK_LEGACY_REMOVE)
352 // We need to expose the enum values at the class level
353 // for backwards compatibility
354 static constexpr CompressionEnum JPEG = CompressionEnum::JPEG;
355 static constexpr CompressionEnum JPEG2000 = CompressionEnum::JPEG2000;
356 static constexpr CompressionEnum JPEGLS = CompressionEnum::JPEGLS;
357 static constexpr CompressionEnum RLE = CompressionEnum::RLE;
358#endif
359
360 itkSetEnumMacro(CompressionType, CompressionEnum);
361 itkGetEnumMacro(CompressionType, CompressionEnum);
362
363 void
364 InternalSetCompressor(const std::string & _compressor) override;
365
366protected:
368 ~GDCMImageIO() override;
369 void
370 PrintSelf(std::ostream & os, Indent indent) const override;
371
372 void
374
375 double m_RescaleSlope{};
376
377 double m_RescaleIntercept{};
378
379 std::string m_UIDPrefix{};
380
381 std::string m_StudyInstanceUID{};
382
383 std::string m_SeriesInstanceUID{};
384
385 std::string m_FrameOfReferenceInstanceUID{};
386
387 bool m_KeepOriginalUID{};
388
389 bool m_LoadPrivateTags{};
390
391 bool m_ReadYBRtoRGB{};
392
393private:
394#if defined(ITKIO_DEPRECATED_GDCM1_API)
395 std::string m_PatientName{};
396
397 std::string m_PatientID{};
398
399 std::string m_PatientDOB{};
400
401 std::string m_StudyID{};
402
403 std::string m_StudyDescription{};
404
405 std::string m_BodyPart{};
406
407 std::string m_NumberOfSeriesInStudy{};
408
409 std::string m_NumberOfStudyRelatedSeries{};
410
411 std::string m_PatientSex{};
412
413 std::string m_PatientAge{};
414
415 std::string m_StudyDate{};
416
417 std::string m_Modality{};
418
419 std::string m_Manufacturer{};
420
421 std::string m_Institution{};
422
423 std::string m_Model{};
424
425 std::string m_ScanOptions{};
426#endif
427
428 unsigned int m_GlobalNumberOfDimensions{};
429
430 CompressionEnum m_CompressionType{};
431
432 bool m_SingleBit{};
433
434 IOComponentEnum m_InternalComponentType{};
435
436 InternalHeader * m_DICOMHeader{};
437};
438
439} // end namespace itk
440
441#endif // itkGDCMImageIO_h
ImageIO class for reading and writing DICOM V3.0 and ACR/NEMA 1&2 uncompressed images....
void WriteImageInformation() override
void InternalReadImageInformation()
void ReadImageInformation() override
~GDCMImageIO() override
void PrintSelf(std::ostream &os, Indent indent) const override
bool CanWriteFile(const char *) override
bool CanReadFile(const char *) override
void InternalSetCompressor(const std::string &_compressor) override
bool GetValueFromTag(const std::string &tag, std::string &value)
static bool GetLabelFromTag(const std::string &tag, std::string &labelId)
void Write(const void *buffer) override
void Read(void *pointer) override
Abstract superclass defines image IO interface.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
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)