DMRITool  v0.1.1-139-g860d86b4
Diffusion MRI Tool
itkMeshFromTensorImageFilter.h
Go to the documentation of this file.
1 
11 #ifndef __itkMeshFromTensorImageFilter_h
12 #define __itkMeshFromTensorImageFilter_h
13 
14 #include <vtkPolyData.h>
15 #include <vtkTensorGlyph.h>
16 #include <vtkDataArray.h>
18 
19 #include "itkDiffusionTensor.h"
20 
21 namespace itk
22 {
23 
33 template < class TInputImage, class TOutputMesh=vtkPolyData>
34 class ITK_EXPORT MeshFromTensorImageFilter
35  : public MeshFromImageImageFilter<TInputImage, TOutputMesh>
36 {
37 public:
41  typedef SmartPointer<Self> Pointer;
42  typedef SmartPointer<const Self> ConstPointer;
43 
45  itkNewMacro(Self);
46 
49 
50  //** Convenient typedefs */
51  typedef TInputImage InputImageType;
52  typedef typename InputImageType::Pointer InputImagePointer;
53  typedef typename InputImageType::ConstPointer InputImageConstPointer;
54  typedef typename InputImageType::IndexType InputImageIndexType;
55  typedef typename InputImageType::SizeType InputImageSizeType;
56  typedef typename InputImageType::SizeValueType InputImageSizeValueType;
57  typedef typename InputImageType::SpacingType InputImageSpacingType;
58  typedef typename InputImageType::PixelType InputImagePixelType;
59  typedef typename InputImageType::RegionType InputImageRegionType;
60  typedef typename InputImageType::PointType InputImagePointType;
61 
62  typedef TOutputMesh OutputMeshPolyDataType;
63 
64  typedef typename Superclass::OutputMeshPointsType OutputMeshPointsType;
65  typedef typename Superclass::OutputMeshCellArrayType OutputMeshCellArrayType;
66  typedef typename Superclass::OutputMeshPolyDataPointer OutputMeshPolyDataPointer;
67  typedef typename Superclass::OutputMeshScalarType OutputMeshScalarType;
68  typedef typename Superclass::OutputMeshRGBType OutputMeshRGBType;
69 
70 
72  typedef typename Superclass::MatrixType MatrixType;
73  typedef typename Superclass::MatrixPointer MatrixPointer;
74  typedef typename Superclass::VectorType VectorType;
75  typedef typename Superclass::VectorPointer VectorPointer;
76  typedef typename Superclass::STDVectorType STDVectorType;
77  typedef typename Superclass::STDVectorPointer STDVectorPointer;
78 
79  typedef Image<double,3> ScalarImageType;
80  typedef typename ScalarImageType::Pointer ScalarImagePointer;
81 
82  typedef vtkTensorGlyph TensorGlyphType;
83 
85  typedef enum
86  {
87  GLYPH_LINE=0, //0
95  GLYPH_SUPERQUADRIC //6
97 
98 
99  typedef enum
100  {
102  COLOR_NONE=0,
110  COLOR_BY_IMAGE
112 
113  itkSetGetMacro(ScalarImage, ScalarImagePointer);
114 
115  itkSetGetMacro(TensorColorScheme, TensorColorSchemeType);
116 
117  itkSetGetMacro(ShapeMode, GlyphShapeType);
118 
119  void SetGlyphResolution(int res);
120  itkGetMacro(GlyphResolution, int);
121 
122  void SetGlyphShape (GlyphShapeType i);
123 
125  void SetGlyphShapeToLine();
126 
128  void SetGlyphShapeToDisk();
129 
131  void SetGlyphShapeToArrow();
132 
134  void SetGlyphShapeToCube();
135 
137  void SetGlyphShapeToCylinder();
138 
140  void SetGlyphShapeToSphere();
141 
143  void SetGlyphShapeToSuperquadric();
144 
145 
146  static double GetScalarCodeFromTensor(const DiffusionTensor<double>& tensor, TensorColorSchemeType colorScheme);
147 
148  // static std::vector<double> GetRGBAFromTensorDirection(const DiffusionTensor<double>& tensor);
149 
150 
151 protected:
153 
155  {};
156 
157  virtual void GenerateData() ITK_OVERRIDE;
158 
159  void PrintSelf(std::ostream& os, Indent indent) const ITK_OVERRIDE
160  {
161  }
162 
163  typename LightObject::Pointer InternalClone() const ITK_OVERRIDE
164  {
165  typename LightObject::Pointer loPtr = Superclass::InternalClone();
166 
167  typename Self::Pointer rval = dynamic_cast<Self *>(loPtr.GetPointer());
168  if(rval.IsNull())
169  {
170  itkExceptionMacro(<< "downcast to type " << this->GetNameOfClass()<< " failed.");
171  }
172 
173  rval->m_TensorColorScheme = m_TensorColorScheme;
174 
175  return loPtr;
176  }
177 
178  void VerifyInputParameters() const ITK_OVERRIDE
179  {
180  if (m_TensorColorScheme==COLOR_BY_IMAGE)
181  {
182  utlSAGlobalException(itk::IsImageEmpty(m_ScalarImage))(this->m_TensorColorScheme).msg("need to set m_ScalarImage for COLOR_BY_IMAGE");
183 
184  InputImageConstPointer inputPtr = this->GetInput();
185  utlSAGlobalException(!itk::VerifyImageSize(inputPtr, m_ScalarImage))
186  .msg("m_ScalarImage and input need to have the same size");
187  }
188  }
189 
190 
191  TensorColorSchemeType m_TensorColorScheme=COLOR_BY_DIRECTION;
192 
193  GlyphShapeType m_ShapeMode = GLYPH_SPHERE;
194 
195  vtkSmartPointer<vtkPolyDataAlgorithm> m_Shape = vtkPolyDataAlgorithm::New();
196 
197  vtkSmartPointer<TensorGlyphType> m_Glyph = TensorGlyphType::New();
198 
199  // vtkSmartPointer<vtkDataArray> m_Scalars = vtkDataArray::New();
200 
202  int m_GlyphResolution=10;
203 
204  ScalarImagePointer m_ScalarImage = ScalarImageType::New();
205 
206 
207 private:
208  MeshFromTensorImageFilter(const Self&); //purposely not implemented
209  void operator=(const Self&);//purposely not implemented
210 
211 };
212 
213 }
214 
215 #if !defined(ITK_MANUAL_INSTANTIATION) && !defined(__itkMeshFromTensorImageFilter_hxx)
217 #endif
218 
219 
220 #endif
void PrintSelf(std::ostream &os, Indent indent) const ITK_OVERRIDE
Superclass::OutputMeshRGBType OutputMeshRGBType
InputImageType::SizeValueType InputImageSizeValueType
InputImageType::PixelType InputImagePixelType
bool IsImageEmpty(const SmartPointer< ImageType > &image)
Definition: utlITK.h:435
Superclass::OutputMeshPolyDataPointer OutputMeshPolyDataPointer
#define utlSAGlobalException(expr)
Definition: utlCoreMacro.h:362
InputImageType::IndexType InputImageIndexType
InputImageType::RegionType InputImageRegionType
#define ITK_OVERRIDE
Definition: utlITKMacro.h:46
Superclass::OutputMeshScalarType OutputMeshScalarType
Superclass::OutputMeshCellArrayType OutputMeshCellArrayType
LightObject::Pointer InternalClone() const ITK_OVERRIDE
Compute mesh from SH coefficients.
Superclass::STDVectorPointer STDVectorPointer
Compute mesh from a tensor image.
#define itkSetGetMacro(name, type)
Definition: utlITKMacro.h:134
void VerifyInputParameters() const ITK_OVERRIDE
InputImageType::PointType InputImagePointType
Superclass::OutputMeshPointsType OutputMeshPointsType
bool VerifyImageSize(const SmartPointer< Image1Type > &image1, const SmartPointer< Image2Type > &image2, const bool isMinimalDimension=false)
Definition: utlITK.h:766
InputImageType::ConstPointer InputImageConstPointer
tensor with some useful functions
MeshFromImageImageFilter< TInputImage, TOutputMesh > Superclass
InputImageType::SpacingType InputImageSpacingType