ANIMA  4.0
animaLogTensorImageFilter.hxx
Go to the documentation of this file.
1 #pragma once
3 
4 #include <itkImageRegionIterator.h>
5 #include <itkImageRegionConstIterator.h>
6 
7 namespace anima
8 {
9 template <class TScalarType, unsigned int NDimensions>
10 void
13 {
14  // Override the method in itkImageSource, so we can set the vector length of
15  // the output itk::VectorImage
16 
17  this->Superclass::GenerateOutputInformation();
18 
19  m_VectorSize = this->GetInput(0)->GetNumberOfComponentsPerPixel();
20  TOutputImage *output = this->GetOutput();
21  output->SetVectorLength(m_VectorSize);
22 }
23 
24 template <class TScalarType, unsigned int NDimensions>
25 void
28 {
29  unsigned int nbInputs = this->GetNumberOfIndexedInputs();
30  if (nbInputs != 1)
31  throw itk::ExceptionObject(__FILE__, __LINE__,"There should be one input... Exiting...",ITK_LOCATION);
32 
33  m_VectorSize = this->GetInput(0)->GetNumberOfComponentsPerPixel();
34  m_TensorDimension = floor((sqrt((double)(8 * m_VectorSize + 1)) - 1) / 2.0);
35 
36  this->GetOutput()->SetNumberOfComponentsPerPixel(m_VectorSize);
37  this->AllocateOutputs();
38 }
39 
40 template <class TScalarType, unsigned int NDimensions>
41 void
43 DynamicThreadedGenerateData (const OutputImageRegionType &outputRegionForThread)
44 {
45  typedef itk::ImageRegionConstIterator< TInputImage > InIteratorType;
46  typedef itk::ImageRegionIterator< TOutputImage > OutRegionIteratorType;
47 
48  InIteratorType inIterator(this->GetInput(), outputRegionForThread);
49  OutRegionIteratorType outIterator(this->GetOutput(), outputRegionForThread);
50 
51  OutputPixelType outValue;
52  vnl_matrix <double> tmpTensor(m_TensorDimension, m_TensorDimension);
53  vnl_matrix <double> tmpLogTensor(m_TensorDimension, m_TensorDimension);
54 
55  LECalculatorPointer leCalculator = LECalculatorType::New();
56 
57  while (!outIterator.IsAtEnd())
58  {
59  outValue = inIterator.Get();
60 
61  if (!isZero(outValue))
62  {
63  anima::GetTensorFromVectorRepresentation(outValue,tmpTensor,m_TensorDimension);
64 
65  leCalculator->GetTensorLogarithm(tmpTensor,tmpLogTensor);
66 
67  anima::GetVectorRepresentation(tmpLogTensor,outValue,m_VectorSize,m_ScaleNonDiagonal);
68  }
69 
70  outIterator.Set(outValue);
71 
72  ++inIterator;
73  ++outIterator;
74  }
75 }
76 
77 } // end of namespace anima
typename LECalculatorType::Pointer LECalculatorPointer
void GetVectorRepresentation(const vnl_matrix< T1 > &tensor, itk::VariableLengthVector< T2 > &vector, unsigned int vecDim=0, bool scale=false)
void BeforeThreadedGenerateData() ITK_OVERRIDE
void GetTensorFromVectorRepresentation(const itk::VariableLengthVector< T1 > &vector, vnl_matrix< T2 > &tensor, unsigned int tensDim=0, bool scale=false)
TOutputImage::PixelType OutputPixelType
Superclass::OutputImageRegionType OutputImageRegionType
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) ITK_OVERRIDE
itk::VectorImage< TScalarType, NDimensions > TOutputImage