File: itkAnisotropicDiffusionVesselEnhancementFunction.txx

package info (click to toggle)
vmtk 1.3%2Bdfsg-2.1%2Bdeb9u1
  • links: PTS, VCS
  • area: non-free
  • in suites: stretch
  • size: 8,932 kB
  • sloc: cpp: 82,947; ansic: 31,817; python: 21,462; perl: 381; makefile: 93; ruby: 41; sh: 19
file content (203 lines) | stat: -rw-r--r-- 6,828 bytes parent folder | download | duplicates (4)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
/*=========================================================================

Program:   Insight Segmentation & Registration Toolkit
Module:    $RCSfile: itkAnisotropicDiffusionVesselEnhancementFunction.txx,v $
Language:  C++
Date:      $Date: 2007/06/20 16:03:23 $
Version:   $Revision: 1.14 $

Copyright (c) Insight Software Consortium. All rights reserved.
See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.

This software is distributed WITHOUT ANY WARRANTY; without even 
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 
PURPOSE.  See the above copyright notices for more information.

=========================================================================*/
#ifndef __itkAnisotropicDiffusionVesselEnhancementFunction_txx
#define __itkAnisotropicDiffusionVesselEnhancementFunction_txx

#include "itkAnisotropicDiffusionVesselEnhancementFunction.h"
#include "vnl/algo/vnl_symmetric_eigensystem.h"

namespace itk {

template< class TImageType >
AnisotropicDiffusionVesselEnhancementFunction< TImageType>
::AnisotropicDiffusionVesselEnhancementFunction()
{
  RadiusType r;

  for( unsigned int i=0 ; i < ImageDimension ; i++ )
    {
    r[i] = 1;
    } 

  this->SetRadius(r);
  
  // Dummy neighborhood.
  NeighborhoodType it;
  it.SetRadius( r );
  
  // Find the center index of the neighborhood.
  m_Center =  it.Size() / 2;

  // Get the stride length for each axis.
  for(unsigned int i = 0; i < ImageDimension; i++)
    {  m_xStride[i] = it.GetStride(i); }
}

template< class TImageType >
typename AnisotropicDiffusionVesselEnhancementFunction< TImageType >
::TimeStepType
AnisotropicDiffusionVesselEnhancementFunction<TImageType>
::ComputeGlobalTimeStep(void *GlobalData) const
{
  /* returns the time step supplied by the user. We don't need
     to use the global data supplied since we are returning a fixed value
  */

  return this->GetTimeStep();
}
 
template< class TImageType >
typename AnisotropicDiffusionVesselEnhancementFunction< TImageType >::PixelType
AnisotropicDiffusionVesselEnhancementFunction< TImageType >
::ComputeUpdate(const NeighborhoodType &it, void *globalData,
                const FloatOffsetType& offset)
{
  double value = 0.0;
  
  return (PixelType) (value);   
}

template< class TImageType >
typename AnisotropicDiffusionVesselEnhancementFunction< TImageType >::PixelType
AnisotropicDiffusionVesselEnhancementFunction< TImageType >
::ComputeUpdate(const NeighborhoodType &it, 
                const DiffusionTensorNeighborhoodType &gt,
                void *globalData,
                const FloatOffsetType& offset)
{
  unsigned int i, j;  
//  const ScalarValueType ZERO = NumericTraits<ScalarValueType>::Zero;
  const ScalarValueType center_value  = it.GetCenterPixel();

  // Global data structure
  GlobalDataStruct *gd = (GlobalDataStruct *)globalData;


  // m_dx -> Intensity first derivative 
  // m_dxy -> Intensity second derivative
  // m_DT_dxy -> Diffusion tensor first derivative

  // Compute the first and 2nd derivative 
  gd->m_GradMagSqr = 1.0e-6;
  for( i = 0; i < ImageDimension; i++)
    {
    const unsigned int positionA = 
      static_cast<unsigned int>( m_Center + m_xStride[i]);
    const unsigned int positionB = 
      static_cast<unsigned int>( m_Center - m_xStride[i]);

    gd->m_dx[i] = 0.5 * (it.GetPixel( positionA ) - 
                     it.GetPixel( positionB )    );

    gd->m_dxy[i][i] = it.GetPixel( positionA )
      + it.GetPixel( positionB ) - 2.0 * center_value;
    
    for( j = i+1; j < ImageDimension; j++ )
      {
      const unsigned int positionAa = static_cast<unsigned int>( 
        m_Center - m_xStride[i] - m_xStride[j] );
      const unsigned int positionBa = static_cast<unsigned int>( 
        m_Center - m_xStride[i] + m_xStride[j] );
      const unsigned int positionCa = static_cast<unsigned int>( 
        m_Center + m_xStride[i] - m_xStride[j] );
      const unsigned int positionDa = static_cast<unsigned int>( 
        m_Center + m_xStride[i] + m_xStride[j] );

      gd->m_dxy[i][j] = gd->m_dxy[j][i] = 0.25 *( it.GetPixel( positionAa )
                                          - it.GetPixel( positionBa )
                                          - it.GetPixel( positionCa )
                                          + it.GetPixel( positionDa )
        );
      }
    }

  // Compute the diffusion tensor matrix first derivatives 
  TensorPixelType center_Tensor_value  = gt.GetCenterPixel();

  for( i = 0; i < ImageDimension; i++)
    {
    const unsigned int positionA = 
      static_cast<unsigned int>( m_Center + m_xStride[i]);
    const unsigned int positionB = 
      static_cast<unsigned int>( m_Center - m_xStride[i]);
    
    TensorPixelType positionA_Tensor_value = gt.GetPixel( positionA );
    TensorPixelType positionB_Tensor_value = gt.GetPixel( positionB );

    for( j = 0; j < ImageDimension; j++)
      { 
      gd->m_DT_dxy[i][j] = 0.5 *  ( positionA_Tensor_value(i,j) - 
                                positionB_Tensor_value(i,j) ); 
      }
    }

  ScalarValueType   pdWrtDiffusion1;

  pdWrtDiffusion1 = gd->m_DT_dxy[0][0] * gd->m_dx[0]  
                    + gd->m_DT_dxy[0][1] * gd->m_dx[1] 
                    + gd->m_DT_dxy[0][2] * gd->m_dx[2];

  ScalarValueType   pdWrtDiffusion2;

  pdWrtDiffusion2 = gd->m_DT_dxy[1][0] * gd->m_dx[0]  
                    + gd->m_DT_dxy[1][1] * gd->m_dx[1] 
                    + gd->m_DT_dxy[1][2] * gd->m_dx[2];

  ScalarValueType  pdWrtDiffusion3;

  pdWrtDiffusion3 = gd->m_DT_dxy[2][0] * gd->m_dx[0]  
                    + gd->m_DT_dxy[2][1] * gd->m_dx[1] 
                    + gd->m_DT_dxy[2][2] * gd->m_dx[2];

  ScalarValueType   pdWrtImageIntensity1;

  pdWrtImageIntensity1 = center_Tensor_value(0,0) *  gd->m_dxy[0][0]  + 
                    center_Tensor_value(0,1) *  gd->m_dxy[0][1] +
                    center_Tensor_value(0,2) *  gd->m_dxy[0][2];
  
  ScalarValueType   pdWrtImageIntensity2;

  pdWrtImageIntensity2 = center_Tensor_value(1,0) *  gd->m_dxy[1][0]  + 
                    center_Tensor_value(1,1) *  gd->m_dxy[1][1] +
                    center_Tensor_value(1,2) *  gd->m_dxy[1][2];
 
  ScalarValueType   pdWrtImageIntensity3;

  pdWrtImageIntensity3 = center_Tensor_value(2,0) *  gd->m_dxy[2][0]  + 
                    center_Tensor_value(2,1) *  gd->m_dxy[2][1] +
                    center_Tensor_value(2,2) *  gd->m_dxy[2][2];
 
  ScalarValueType   total;

  total = pdWrtDiffusion1 + pdWrtDiffusion2 + pdWrtDiffusion3 +
         pdWrtImageIntensity1 + pdWrtImageIntensity2 + pdWrtImageIntensity3;

  return ( PixelType ) ( total );
} 

template <class TImageType>
void
AnisotropicDiffusionVesselEnhancementFunction<TImageType>::
PrintSelf(std::ostream& os, Indent indent) const
{
  Superclass::PrintSelf(os, indent);
}


} // end namespace itk

#endif