ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkImageIOBase.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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 #ifndef itkImageIOBase_h
19 #define itkImageIOBase_h
20 #include "ITKIOImageBaseExport.h"
21 
22 #include "itkIOConfigure.h"
23 
24 #include "itkLightProcessObject.h"
25 #include "itkIndent.h"
26 #include "itkImageIORegion.h"
27 #include "itkRGBPixel.h"
28 #include "itkRGBAPixel.h"
29 #include "itkCovariantVector.h"
31 #include "itkDiffusionTensor3D.h"
33 
34 #include "vnl/vnl_vector.h"
35 #include "vcl_compiler.h"
36 
37 #include <fstream>
38 #include <string>
39 
40 namespace itk
41 {
42 // Forward reference for VariableLengthVector
43 template <typename TValue> class VariableLengthVector;
44 
73 class ITKIOImageBase_EXPORT ImageIOBase:public LightProcessObject
74 {
75 public:
77  typedef ImageIOBase Self;
80 
82  itkTypeMacro(ImageIOBase, Superclass);
83 
85  itkSetStringMacro(FileName);
86  itkGetStringMacro(FileName);
88 
92 
98  class UnknownType {};
99 
103  typedef enum { UNKNOWNPIXELTYPE, SCALAR, RGB, RGBA, OFFSET, VECTOR,
104  POINT, COVARIANTVECTOR, SYMMETRICSECONDRANKTENSOR,
105  DIFFUSIONTENSOR3D, COMPLEX, FIXEDARRAY, MATRIX } IOPixelType;
106 
111  typedef enum { UNKNOWNCOMPONENTTYPE, UCHAR, CHAR, USHORT, SHORT, UINT, INT,
112  ULONG, LONG, FLOAT, DOUBLE } IOComponentType;
113 
117  void SetNumberOfDimensions(unsigned int);
118 
119  itkGetConstMacro(NumberOfDimensions, unsigned int);
120 
124  virtual void SetDimensions(unsigned int i, unsigned int dim);
126 
127  virtual itk::SizeValueType GetDimensions(unsigned int i) const
128  { return m_Dimensions[i]; }
129 
132  virtual void SetOrigin(unsigned int i, double origin);
133 
134  virtual double GetOrigin(unsigned int i) const
135  {
136  return m_Origin[i];
137  }
138 
141  virtual void SetSpacing(unsigned int i, double spacing);
142 
143  virtual double GetSpacing(unsigned int i) const
144  {
145  return m_Spacing[i];
146  }
147 
150  virtual void SetDirection(unsigned int i, const std::vector< double > & direction);
151 
152  virtual void SetDirection(unsigned int i, const vnl_vector< double > & direction);
153 
154  virtual std::vector< double > GetDirection(unsigned int i) const
155  {
156  return m_Direction[i];
157  }
158 
161  virtual std::vector< double > GetDefaultDirection(unsigned int i) const;
162 
169  itkSetMacro(IORegion, ImageIORegion);
170  itkGetConstReferenceMacro(IORegion, ImageIORegion);
172 
178  itkSetEnumMacro(PixelType, IOPixelType);
179  itkGetEnumMacro(PixelType, IOPixelType);
181 
184  itkSetEnumMacro(ComponentType, IOComponentType);
185  itkGetEnumMacro(ComponentType, IOComponentType);
186 
192  virtual const std::type_info & GetComponentTypeInfo() const;
193 
198  itkSetMacro(NumberOfComponents, unsigned int);
199  itkGetConstReferenceMacro(NumberOfComponents, unsigned int);
201 
203  itkSetMacro(UseCompression, bool);
204  itkGetConstMacro(UseCompression, bool);
205  itkBooleanMacro(UseCompression);
207 
209  itkSetMacro(UseStreamedReading, bool);
210  itkGetConstMacro(UseStreamedReading, bool);
211  itkBooleanMacro(UseStreamedReading);
213 
215  itkSetMacro(UseStreamedWriting, bool);
216  itkGetConstMacro(UseStreamedWriting, bool);
217  itkBooleanMacro(UseStreamedWriting);
219 
222  static std::string GetComponentTypeAsString(IOComponentType);
223 
225  static IOComponentType GetComponentTypeFromString(const std::string &typeString);
226 
229  static std::string GetPixelTypeAsString(IOPixelType);
230 
232  static IOPixelType GetPixelTypeFromString(const std::string &pixelString);
233 
236  typedef enum { ASCII, Binary, TypeNotApplicable } FileType;
237 
240  typedef enum { BigEndian, LittleEndian, OrderNotApplicable } ByteOrder;
241 
244  itkSetEnumMacro(FileType, FileType);
245  itkGetEnumMacro(FileType, FileType);
247  {
248  this->SetFileType(ASCII);
249  }
251 
253  {
254  this->SetFileType(Binary);
255  }
256 
268  itkSetEnumMacro(ByteOrder, ByteOrder);
269  itkGetEnumMacro(ByteOrder, ByteOrder);
271  {
272  this->SetByteOrder(BigEndian);
273  }
275 
277  {
278  this->SetByteOrder(LittleEndian);
279  }
280 
283  std::string GetFileTypeAsString(FileType) const;
284 
287  std::string GetByteOrderAsString(ByteOrder) const;
288 
291 
295 
302  virtual SizeType GetPixelStride() const;
303 
305  SizeType GetImageSizeInPixels() const;
306 
308  SizeType GetImageSizeInBytes() const;
309 
312  SizeType GetImageSizeInComponents() const;
313 
318  virtual unsigned int GetComponentSize() const;
319 
320  /*-------- This part of the interfaces deals with reading data ----- */
321 
324  virtual bool CanReadFile(const char *) = 0;
325 
330  virtual bool CanStreamRead()
331  {
332  return false;
333  }
334 
337  virtual void ReadImageInformation() = 0;
338 
340  virtual void Read(void *buffer) = 0;
341 
342  /*-------- This part of the interfaces deals with writing data ----- */
343 
346  virtual bool CanWriteFile(const char *) = 0;
347 
354  virtual bool CanStreamWrite()
355  {
356  return false;
357  }
358 
361  virtual void WriteImageInformation() = 0;
362 
366  virtual void Write(const void *buffer) = 0;
367 
368  /* --- Support reading and writing data as a series of files. --- */
369 
375  virtual bool SupportsDimension(unsigned long dim)
376  {
377  return ( dim == 2 );
378  }
379 
391  virtual ImageIORegion
392  GenerateStreamableReadRegionFromRequestedRegion(const ImageIORegion & requested) const;
393 
408  virtual unsigned int GetActualNumberOfSplitsForWriting(unsigned int numberOfRequestedSplits,
409  const ImageIORegion & pasteRegion,
410  const ImageIORegion & largestPossibleRegion);
411 
418  virtual ImageIORegion GetSplitRegionForWriting(unsigned int ithPiece,
419  unsigned int numberOfActualSplits,
420  const ImageIORegion & pasteRegion,
421  const ImageIORegion & largestPossibleRegion);
422 
424  typedef std::vector< std::string > ArrayOfExtensionsType;
425 
430  const ArrayOfExtensionsType & GetSupportedReadExtensions() const;
431 
436  const ArrayOfExtensionsType & GetSupportedWriteExtensions() const;
437 
438  template <typename TPixel>
439  void SetTypeInfo(const TPixel *);
440 
442  template <typename TPixel>
444  {
445  static const IOComponentType CType =
446  UNKNOWNCOMPONENTTYPE;
447  };
448  template <typename TPixel>
449  void SetPixelTypeInfo(const TPixel *)
450  {
451  this->SetNumberOfComponents(1);
452  this->SetPixelType(SCALAR);
453  this->SetComponentType(MapPixelType<TPixel>::CType);
454  }
455  template <typename TPixel>
457  {
458  this->SetNumberOfComponents(3);
459  this->SetPixelType(RGB);
460  this->SetComponentType(MapPixelType<TPixel>::CType);
461  }
462  template <typename TPixel>
464  {
465  this->SetNumberOfComponents(4);
466  this->SetPixelType(RGBA);
467  this->SetComponentType(MapPixelType<TPixel>::CType);
468  }
469  template <typename TPixel, unsigned VLength>
471  {
472  this->SetNumberOfComponents(VLength);
473  this->SetPixelType(VECTOR);
474  this->SetComponentType(MapPixelType<TPixel>::CType);
475  }
476  template <typename TPixel>
478  {
479  this->SetNumberOfComponents(1);
480  this->SetPixelType(VECTOR);
481  this->SetComponentType(MapPixelType<TPixel>::CType);
482  }
483  template <typename TPixel, unsigned VLength>
485  {
486  this->SetNumberOfComponents(VLength);
487  this->SetPixelType(COVARIANTVECTOR);
488  this->SetComponentType(MapPixelType<TPixel>::CType);
489  }
490  template <typename TPixel,unsigned VLength>
492  {
493  this->SetNumberOfComponents(VLength);
494  this->SetPixelType(COVARIANTVECTOR);
495  this->SetComponentType(MapPixelType<TPixel>::CType);
496  }
498 
499  template <typename TPixel, unsigned VLength>
501  {
502  this->SetNumberOfComponents(VLength * (VLength + 1) / 2 );
503  this->SetPixelType(SYMMETRICSECONDRANKTENSOR);
504  this->SetComponentType(MapPixelType<TPixel>::CType);
505  }
506 
507  template <typename TPixel>
509  {
510  this->SetNumberOfComponents(6);
511  this->SetPixelType(DIFFUSIONTENSOR3D);
512  this->SetComponentType(MapPixelType<TPixel>::CType);
513  }
514 
515  template <typename TPixel, unsigned VLength>
517  {
518  this->SetNumberOfComponents(VLength * VLength);
519  this->SetPixelType(MATRIX);
520  this->SetComponentType(MapPixelType<TPixel>::CType);
521  }
522 
523  template <typename TPixel>
524  void SetPixelTypeInfo(const std::complex< TPixel > *)
525  {
526  this->SetNumberOfComponents(2);
527  this->SetPixelType(COMPLEX);
528  this->SetComponentType(MapPixelType<TPixel>::CType);
529  }
530 
531  template <unsigned VLength>
533  {
534  this->SetNumberOfComponents(VLength);
535  this->SetPixelType(ImageIOBase::OFFSET);
536  this->SetComponentType(ImageIOBase::LONG);
537  }
538 
539 protected:
540  ImageIOBase();
541  ~ImageIOBase();
542  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
543 
544  virtual const ImageRegionSplitterBase* GetImageRegionSplitter() const;
545 
548 
552 
555 
557 
560 
562  std::string m_FileName;
563 
566  unsigned int m_NumberOfComponents;
567 
569  unsigned int m_NumberOfDimensions;
570 
573 
576 
579 
583 
585  std::vector< SizeValueType > m_Dimensions;
586 
589  std::vector< double > m_Spacing;
590 
592  std::vector< double > m_Origin;
593 
595  std::vector< std::vector< double > > m_Direction;
596 
599  std::vector< SizeType > m_Strides;
600 
602  virtual void Reset(const bool freeDynamic = true);
603 
605  void Resize(const unsigned int numDimensions,
606  const unsigned int *dimensions);
607 
610  virtual unsigned int GetPixelSize() const;
611 
618  void ComputeStrides();
619 
622  SizeType GetComponentStride() const;
623 
626  SizeType GetRowStride() const;
627 
630  SizeType GetSliceStride() const;
631 
643  virtual void OpenFileForReading(std::ifstream & inputStream, const std::string & filename,
644  bool ascii = false);
645 
661  virtual void OpenFileForWriting(std::ofstream & outputStream, const std::string & filename,
662  bool truncate = true, bool ascii = false);
663 
665  virtual void WriteBufferAsASCII(std::ostream & os, const void *buffer,
666  IOComponentType ctype,
667  SizeType numberOfBytesToWrite);
668 
670  virtual void ReadBufferAsASCII(std::istream & os, void *buffer,
671  IOComponentType ctype,
672  SizeType numberOfBytesToBeRead);
673 
675  bool ReadBufferAsBinary(std::istream & os, void *buffer, SizeType numberOfBytesToBeRead);
676 
678  void AddSupportedReadExtension(const char *extension);
679 
681  void AddSupportedWriteExtension(const char *extension);
682 
685  virtual unsigned int GetActualNumberOfSplitsForWritingCanStreamWrite(unsigned int numberOfRequestedSplits,
686  const ImageIORegion & pasteRegion) const;
687 
690  virtual ImageIORegion GetSplitRegionForWritingCanStreamWrite(unsigned int ithPiece,
691  unsigned int numberOfActualSplits,
692  const ImageIORegion & pasteRegion) const;
693 
694 private:
695  ImageIOBase(const Self &) ITK_DELETE_FUNCTION;
696  void operator=(const Self &) ITK_DELETE_FUNCTION;
697 
698  ArrayOfExtensionsType m_SupportedReadExtensions;
699  ArrayOfExtensionsType m_SupportedWriteExtensions;
700 };
701 
702 #define IMAGEIOBASE_TYPEMAP(type,ctype) \
703  template <> struct ImageIOBase::MapPixelType<type> \
704  { \
705  static const IOComponentType CType = ctype; \
706  }
707 
708 // the following typemaps are not platform independent
709 #if VCL_CHAR_IS_SIGNED
710 IMAGEIOBASE_TYPEMAP(signed char, CHAR);
711 #endif // VCL_CHAR_IS_SIGNED
712 IMAGEIOBASE_TYPEMAP(char, CHAR);
713 IMAGEIOBASE_TYPEMAP(unsigned char, UCHAR);
714 IMAGEIOBASE_TYPEMAP(short, SHORT);
715 IMAGEIOBASE_TYPEMAP(unsigned short, USHORT);
716 IMAGEIOBASE_TYPEMAP(int, INT);
717 IMAGEIOBASE_TYPEMAP(unsigned int, UINT);
718 IMAGEIOBASE_TYPEMAP(long, LONG);
719 IMAGEIOBASE_TYPEMAP(unsigned long, ULONG);
720 IMAGEIOBASE_TYPEMAP(float, FLOAT);
721 IMAGEIOBASE_TYPEMAP(double, DOUBLE);
722 #undef IMAGIOBASE_TYPEMAP
723 
724 
725 } // end namespace itk
726 
727 #endif // itkImageIOBase_h
void SetPixelTypeInfo(const Offset< VLength > *)
A templated class holding a M x N size Matrix.
Definition: itkMatrix.h:47
unsigned int m_NumberOfComponents
Light weight base class for most itk classes.
virtual itk::SizeValueType GetDimensions(unsigned int i) const
ImageIOBase Self
Represent the offset between two n-dimensional indexes in a n-dimensional image.
Definition: itkOffset.h:55
virtual bool SupportsDimension(unsigned long dim)
int64_t intmax_t
Definition: itkIntTypes.h:114
void SetPixelTypeInfo(const DiffusionTensor3D< TPixel > *)
An ImageIORegion represents a structured region of data.
void SetByteOrderToLittleEndian()
Abstract superclass defines image IO interface.
signed long OffsetValueType
Definition: itkIntTypes.h:154
Represent Red, Green, Blue and Alpha components for color images.
Definition: itkRGBAPixel.h:59
Represent a symmetric tensor of second rank.
ByteOrder m_ByteOrder
signed long IndexValueType
Definition: itkIntTypes.h:150
void SetPixelTypeInfo(const Matrix< TPixel, VLength, VLength > *)
ImageIORegion m_IORegion
std::vector< double > m_Origin
void SetPixelTypeInfo(const VariableLengthVector< TPixel > *)
std::vector< double > m_Spacing
LightProcessObject Superclass
virtual double GetSpacing(unsigned int i) const
std::vector< std::string > ArrayOfExtensionsType
unsigned long SizeValueType
Definition: itkIntTypes.h:143
virtual std::vector< double > GetDirection(unsigned int i) const
IOPixelType m_PixelType
void SetFileTypeToBinary()
unsigned int m_NumberOfDimensions
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
void SetPixelTypeInfo(const Vector< TPixel, VLength > *)
void SetPixelTypeInfo(const TPixel *)
::itk::intmax_t SizeType
virtual double GetOrigin(unsigned int i) const
Represents an array whose length can be defined at run-time.
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
std::vector< std::vector< double > > m_Direction
Divide an image region into several pieces.
void SetPixelTypeInfo(const SymmetricSecondRankTensor< TPixel, VLength > *)
void SetPixelTypeInfo(const std::complex< TPixel > *)
virtual bool CanStreamRead()
::itk::SizeValueType SizeValueType
IOComponentType m_ComponentType
void SetPixelTypeInfo(const CovariantVector< TPixel, VLength > *)
::itk::OffsetValueType BufferSizeType
void SetPixelTypeInfo(const RGBPixel< TPixel > *)
::itk::IndexValueType IndexValueType
std::vector< SizeType > m_Strides
LightProcessObject is the base class for all process objects (source, filters, mappers) in the Insigh...
std::string m_FileName
Control indentation during Print() invocation.
Definition: itkIndent.h:49
IMAGEIOBASE_TYPEMAP(char, CHAR)
SmartPointer< Self > Pointer
void SetByteOrderToBigEndian()
Base class for most ITK classes.
Definition: itkObject.h:57
Represent a diffusion tensor as used in DTI images.
A templated class holding a n-Dimensional covariant vector.
void SetPixelTypeInfo(const RGBAPixel< TPixel > *)
std::vector< SizeValueType > m_Dimensions
void SetPixelTypeInfo(const FixedArray< TPixel, VLength > *)
virtual bool CanStreamWrite()