Shapeworks Studio  2.1
Shape analysis software suite
itkPSMImageDomainWithCurvature.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 __itkPSMImageDomainWithCurvature_h
19 #define __itkPSMImageDomainWithCurvature_h
20 
21 #include "itkPSMImageDomainWithHessians.h"
22 #include "itkImageRegionIteratorWithIndex.h"
23 #include "itkImageRegionIterator.h"
24 #include "itkDiscreteGaussianImageFilter.h"
25 
26 namespace itk
27 {
39 template <class T, unsigned int VDimension>
41  : public PSMImageDomainWithHessians<T, VDimension>
42 {
43 public:
47  typedef SmartPointer<Self> Pointer;
48  typedef SmartPointer<const Self> ConstPointer;
49  typedef WeakPointer<const Self> ConstWeakPointer;
50 
51  typedef typename Superclass::PointType PointType;
52  typedef typename Superclass::ImageType ImageType;
53  typedef typename Superclass::ScalarInterpolatorType ScalarInterpolatorType;
54  typedef typename Superclass::VnlMatrixType VnlMatrixType;
55 
57  itkNewMacro(Self);
58 
61 
63  itkStaticConstMacro(Dimension, unsigned int, VDimension);
64 
67  void SetImage(ImageType *I)
68  {
69  // Computes partial derivatives in parent class
70  Superclass::SetImage(I);
71 
72  typedef itk::DiscreteGaussianImageFilter<ImageType, ImageType> DiscreteGaussianImageFilterType;
73  typename DiscreteGaussianImageFilterType::Pointer f = DiscreteGaussianImageFilterType::New();
74 
75  double sig = this->GetImage()->GetSpacing()[0] * 0.5;
76 
77  f->SetVariance(sig);
78  f->SetInput(this->GetImage());
79  f->SetUseImageSpacingOn();
80  f->Update();
81 
82  // NOTE: Running the image through a filter seems to be the only way
83  // to get all of the image information correct! Set the variance to
84  // positive value to smooth the curvature calculations.
85  m_CurvatureImage = f->GetOutput();
86 
87  // Loop through the image and compute all curvature values
88  itk::ImageRegionIteratorWithIndex<ImageType> it(f->GetOutput(),
89  f->GetOutput()->GetRequestedRegion());
90  itk::ImageRegionIterator<ImageType> oit(m_CurvatureImage, f->GetOutput()->GetRequestedRegion());
91 
92  for (; !it.IsAtEnd(); ++it, ++oit)
93  {
94  // Only compute in a narrow band.
95  if (it.Get() < 4.0 && it.Get() > -4.0)
96  {
97  PointType pos;
98  this->GetImage()->TransformIndexToPhysicalPoint(it.GetIndex(), pos);
99  oit.Set(this->MeanCurvature(pos));
100  }
101  else oit.Set(1.0e-6);
102  }
103 
104  // Release the memory in the parent hessian images.
105  this->DeletePartialDerivativeImages();
106 
107  m_CurvatureInterpolator->SetInputImage(m_CurvatureImage);
108  } // end setimage
109 
110  double GetCurvature(const PointType &pos) const
111  {
112  return m_CurvatureInterpolator->Evaluate(pos);
113  }
114 
115  typename ImageType::Pointer *GetCurvatureImage()
116  { return m_CurvatureImage; }
117 
118 protected:
120  {
121  m_CurvatureInterpolator = ScalarInterpolatorType::New();
122  m_CurvatureImage = ImageType::New();
123  }
124 
125  void PrintSelf(std::ostream& os, Indent indent) const
126  {
127  Superclass::PrintSelf(os, indent);
128  }
129  virtual ~PSMImageDomainWithCurvature() {};
130 
131  double MeanCurvature(const PointType& pos)
132  {
133  // See Kindlmann paper "Curvature-Based Transfer Functions for Direct Volume
134  // Rendering..." for detailss
135 
136  // Get the normal vector associated with this position.
137  //VnlVectorType posnormal = this->SampleNormalVnl(pos, 1.0e-10);
138  typename Superclass::VnlVectorType posnormal = this->SampleNormalVnl(pos, 1.0e-6);
139 
140  // Sample the Hessian for this point and compute gradient of the normal.
141  typename Superclass::VnlMatrixType I;
142  I.set_identity();
143 
144  typename Superclass::VnlMatrixType H = this->SampleHessianVnl(pos);
145  typename Superclass::VnlMatrixType P = I - outer_product(posnormal, posnormal);
146  typename Superclass::VnlMatrixType G = P.transpose() * H * P;
147 
148  // Compute the principal curvatures k1 and k2
149  double frobnorm = sqrt(fabs(this->vnl_trace(G * G.transpose())) + 1.0e-10);
150  double trace = this->vnl_trace(G);
151  // std::cout << "trace = " << trace << std::endl;
152  // std::cout << "G = " << G << std::endl;
153  // std::cout << "frobnorm = " << frobnorm << std::endl;
154  double diff1 = frobnorm * frobnorm *2.0;
155  double diff2 = trace * trace;
156  double frobnorm2;
157 
158  if (frobnorm <= trace) frobnorm2 = 0.0;
159  else frobnorm2 = sqrt((diff1 - diff2 < 0.0) ? 0.0 : diff1 - diff2);
160  double k1 = (trace + frobnorm2) * 0.5;
161  double k2 = (trace - frobnorm2) * 0.5;
162  // std::cout << "k1 = " << k1 << std::endl;
163  // std::cout << "k2= " << k2 << std::endl;
164  // Compute the mean curvature.
165  // double mc = fabs((k1 + k2) / 2.0);
166  return sqrt(k1 * k1 + k2 * k2);
167  }
168 
169  inline double vnl_trace(const VnlMatrixType &m) const
170  {
171  double sum = 0.0;
172  for (unsigned int i = 0; i < VDimension; i++)
173  {
174  sum += m[i][i];
175  }
176  return sum;
177  }
178 
179 private:
180  PSMImageDomainWithCurvature(const Self&); //purposely not implemented
181  void operator=(const Self&); //purposely not implemented
182 
183  // Curvature values are stored in an image
184  typename ImageType::Pointer m_CurvatureImage;
185  typename ScalarInterpolatorType::Pointer m_CurvatureInterpolator;
186 };
187 
188 } // end namespace itk
189 
190 #endif
Image< T, VDimension > ImageType