DMRITool  v0.1.1-139-g860d86b4
Diffusion MRI Tool
itkSpatiallyDenseSparseVectorImageFileWriter.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Spatially Dense Sparse Vector Image File Writer
4 
5  Copyright (c) Pew-Thian Yap. All rights reserved.
6  See http://www.unc.edu/~ptyap/ for details.
7 
8  This software is distributed WITHOUT ANY WARRANTY; without even
9  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
10  PURPOSE. See the above copyright notices for more information.
11 
12  =========================================================================*/
13 
14 #ifndef __itkSpatiallyDenseSparseVectorImageFileWriter_h
15 #define __itkSpatiallyDenseSparseVectorImageFileWriter_h
16 
17 #include "itkProcessObject.h"
18 #include "itkExceptionObject.h"
20 #include "itkImage.h"
21 #include "itkImageFileWriter.h"
22 #include "itkMacro.h"
23 
24 namespace itk
25 {
26 
34 template <class TInputImage>
35 class ITKIOImageBase_HIDDEN SpatiallyDenseSparseVectorImageFileWriter : public ProcessObject
36 {
37 public:
40  typedef ProcessObject Superclass;
41  typedef SmartPointer<Self> Pointer;
42  typedef SmartPointer<const Self> ConstPointer;
43 
45  itkNewMacro(Self);
46 
48  itkTypeMacro(SpatiallyDenseSparseVectorImageFileWriter,ProcessObject);
49 
51  typedef TInputImage InputImageType;
52  typedef typename InputImageType::Pointer InputImagePointer;
53  typedef typename InputImageType::RegionType InputImageRegionType;
54  typedef typename InputImageType::SizeType InputImageSizeType;
55  typedef typename InputImageType::IndexType InputImageIndexType;
56  typedef typename InputImageType::SpacingType InputImageSpacingType;
57  typedef typename InputImageType::PointType InputImagePointType;
58  typedef typename InputImageType::DirectionType InputImageDirectionType;
59  typedef typename InputImageType::InternalPixelType InputImageInternalPixelType;
60  typedef typename InputImageType::ValueType InputImageValueType;
61  typedef typename InputImageType::KeyType InputImageKeyType;
62  typedef typename InputImageType::PixelContainer InputImagePixelContainerType;
63  typedef typename InputImageType::PixelContainerConstPointer InputImagePixelContainerConstPointerType;
64  typedef typename InputImageType::PixelMapType InputImagePixelPixelMapType;
65  typedef typename InputImageType::PixelMapConstIterator InputImagePixelMapConstIteratorType;
66 
67  typedef Image<InputImageKeyType, 1> KeyImageType;
68  typedef typename KeyImageType::Pointer KeyImagePointer;
69  typedef Image<InputImageValueType, 1> ValueImageType;
70  typedef typename ValueImageType::Pointer ValueImagePointer;
71 
72  typedef ImageFileWriter<KeyImageType> KeyImageFileWriterType;
73  typedef typename KeyImageFileWriterType::Pointer KeyImageFileWriterPointer;
74  typedef ImageFileWriter<ValueImageType> ValueImageFileWriterType;
75  typedef typename ValueImageFileWriterType::Pointer ValueImageFileWriterPointer;
76 
78  using Superclass::SetInput;
79  void SetInput(const InputImageType *input);
80  const InputImageType * GetInput(void);
81  const InputImageType * GetInput(unsigned int idx);
82 
84  itkSetStringMacro(FileName);
85  itkGetStringMacro(FileName);
86 
94  virtual void Write(void);
95 
98  virtual void Update() ITK_OVERRIDE
99  {
100  this->Write();
101  }
102 
104  itkSetMacro(UseCompression,bool);
105  itkGetConstReferenceMacro(UseCompression,bool);
106  itkBooleanMacro(UseCompression);
107 
114 // itkSetMacro(UseInputMetaDataDictionary,bool);
115 // itkGetConstReferenceMacro(UseInputMetaDataDictionary,bool);
116 // itkBooleanMacro(UseInputMetaDataDictionary);
117 
118 
119 protected:
122  void PrintSelf(std::ostream& os, Indent indent) const ITK_OVERRIDE;
123 
124  KeyImagePointer m_KeyImage;
125  ValueImagePointer m_ValueImage;
126 // LengthImagePointer m_LengthImage;
127 
128  KeyImageFileWriterPointer m_KeyImageFileWriter;
129  ValueImageFileWriterPointer m_ValueImageFileWriter;
130 // LengthImageFileWriterPointer m_ValueImageFileWriter;
131 
133  void GenerateData(void) ITK_OVERRIDE;
134 
135 private:
136  SpatiallyDenseSparseVectorImageFileWriter(const Self&); //purposely not implemented
137  void operator=(const Self&); //purposely not implemented
138 
139  std::string m_FileName;
140 
141  bool m_UseCompression;
142 // bool m_UseInputMetaDataDictionary; // whether to use the
143  // MetaDataDictionary from the
144  // input or not.
145 };
146 
147 
148 } // end namespace itk
149 
150 #ifndef ITK_MANUAL_INSTANTIATION
152 #endif
153 
154 #endif // __itkSpatiallyDenseSparseVectorImageFileWriter_h
STL namespace.
Writes sparse vector image data to key and value files.
#define ITK_OVERRIDE
Definition: utlITKMacro.h:46
InputImageType::PixelContainerConstPointer InputImagePixelContainerConstPointerType