Shapeworks Studio  2.1
Shape analysis software suite
itkPSMImageDomainWithHessians.h
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 __itkPSMImageDomainWithHessians_h
19 #define __itkPSMImageDomainWithHessians_h
20 
21 #include "itkImage.h"
22 #include "itkPSMImageDomainWithGradients.h"
23 #include "itkLinearInterpolateImageFunction.h"
24 #include "itkGradientImageFilter.h"
25 #include "itkFixedArray.h"
26 #include "itkImageDuplicator.h"
27 #include "itkDiscreteGaussianImageFilter.h"
28 #include "itkDerivativeImageFilter.h"
29 #include "vnl/vnl_matrix_fixed.h"
30 
31 namespace itk
32 {
44 template <class T, unsigned int VDimension>
46  : public PSMImageDomainWithGradients<T, VDimension>
47 {
48 public:
52  typedef SmartPointer<Self> Pointer;
53  typedef SmartPointer<const Self> ConstPointer;
54  typedef WeakPointer<const Self> ConstWeakPointer;
55 
57  typedef typename Superclass::PointType PointType;
58 
59  typedef typename Superclass::ImageType ImageType;
60  typedef typename Superclass::ScalarInterpolatorType ScalarInterpolatorType;
61  typedef vnl_matrix_fixed<T, VDimension, VDimension> VnlMatrixType;
62 
64  itkNewMacro(Self);
65 
68 
70  itkStaticConstMacro(Dimension, unsigned int, VDimension);
71 
74  void SetImage(ImageType *I)
75  {
76  Superclass::SetImage(I);
77 
78  typename DiscreteGaussianImageFilter<ImageType, ImageType>::Pointer
79  gaussian = DiscreteGaussianImageFilter<ImageType, ImageType>::New();
80  gaussian->SetVariance(m_Sigma * m_Sigma);
81  gaussian->SetInput(this->GetImage());
82  gaussian->SetUseImageSpacingOn();
83  gaussian->Update();
84 
85  // Compute the second derivatives and set up the interpolators
86  for (unsigned int i = 0; i < VDimension; i++)
87  {
88  typename DerivativeImageFilter<ImageType, ImageType>::Pointer
89  deriv = DerivativeImageFilter<ImageType, ImageType>::New();
90  deriv->SetInput(gaussian->GetOutput());
91  deriv->SetDirection(i);
92  deriv->SetOrder(2);
93  deriv->SetUseImageSpacingOn();
94  deriv->Update();
95 
96  m_PartialDerivatives[i] = deriv->GetOutput();
97 
98  m_Interpolators[i] = ScalarInterpolatorType::New();
99  m_Interpolators[i]->SetInputImage(m_PartialDerivatives[i]);
100  }
101 
102  // Compute the cross derivatives and set up the interpolators
103  unsigned int k = VDimension;
104  for (unsigned int i = 0; i < VDimension; i++)
105  {
106  for (unsigned int j = i+1; j < VDimension; j++, k++)
107  {
108  typename DerivativeImageFilter<ImageType, ImageType>::Pointer
109  deriv1 = DerivativeImageFilter<ImageType, ImageType>::New();
110  deriv1->SetInput(gaussian->GetOutput());
111  deriv1->SetDirection(i);
112  deriv1->SetUseImageSpacingOn();
113  deriv1->SetOrder(1);
114  deriv1->Update();
115 
116  typename DerivativeImageFilter<ImageType, ImageType>::Pointer
117  deriv2 = DerivativeImageFilter<ImageType, ImageType>::New();
118  deriv2->SetInput(deriv1->GetOutput());
119  deriv2->SetDirection(j);
120  deriv2->SetUseImageSpacingOn();
121  deriv2->SetOrder(1);
122 
123  deriv2->Update();
124 
125  m_PartialDerivatives[k] = deriv2->GetOutput();
126  m_Interpolators[k] = ScalarInterpolatorType::New();
127  m_Interpolators[k]->SetInputImage(m_PartialDerivatives[k]);
128  }
129  }
130  } // end setimage
131 
135  inline VnlMatrixType SampleHessianVnl(const PointType &p) const
136  {
137  VnlMatrixType ans;
138  for (unsigned int i = 0; i < VDimension; i++)
139  { ans[i][i] = m_Interpolators[i]->Evaluate(p); }
140 
141  // Cross derivatives
142  unsigned int k = VDimension;
143  for (unsigned int i =0; i < VDimension; i++)
144  {
145  for (unsigned int j = i+1; j < VDimension; j++, k++)
146  {
147  ans[i][j] = ans[j][i] = m_Interpolators[k]->Evaluate(p);
148  }
149  }
150  return ans;
151  }
152 
157  itkSetMacro(Sigma, double);
158  itkGetMacro(Sigma, double);
159 
161  typename ScalarInterpolatorType::Pointer *GetInterpolators()
162  { return m_Interpolators; }
163  typename ImageType::Pointer *GetPartialDerivatives()
164  { return m_PartialDerivatives; }
165 
166 protected:
167  PSMImageDomainWithHessians() : m_Sigma(0.0)
168  { }
169 
170  void PrintSelf(std::ostream& os, Indent indent) const
171  {
172  Superclass::PrintSelf(os, indent);
173  }
174  virtual ~PSMImageDomainWithHessians() {};
175 
176  void DeletePartialDerivativeImages()
177  {
178  for (unsigned int i = 0; i < VDimension + ((VDimension * VDimension) - VDimension) / 2; i++)
179  {
180  m_PartialDerivatives[i]=0;
181  m_Interpolators[i]=0;
182  }
183  }
184 
185 private:
186  double m_Sigma;
187 
188  PSMImageDomainWithHessians(const Self&); //purposely not implemented
189  void operator=(const Self&); //purposely not implemented
190 
191  // Partials are stored:
192  // 0: dxx 3: dxy 4: dxz
193  // 1: dyy 5: dyz
194  // 2: dzz
195  //
196  typename ImageType::Pointer m_PartialDerivatives[ VDimension + ((VDimension * VDimension) - VDimension) / 2];
197 
198  typename ScalarInterpolatorType::Pointer m_Interpolators[VDimension + ((VDimension * VDimension) - VDimension) / 2];
199 };
200 
201 } // end namespace itk
202 
203 #endif
Image< T, VDimension > ImageType
VnlMatrixType SampleHessianVnl(const PointType &p) const
ScalarInterpolatorType::Pointer * GetInterpolators()