Medical Imaging Interaction Toolkit  2016.11.0
Medical Imaging Interaction Toolkit
mitkCovarianceMatrixCalculator.cpp
Go to the documentation of this file.
1 /*===================================================================
2 
3 The Medical Imaging Interaction Toolkit (MITK)
4 
5 Copyright (c) German Cancer Research Center,
6 Division of Medical and Biological Informatics.
7 All rights reserved.
8 
9 This software is distributed WITHOUT ANY WARRANTY; without
10 even the implied warranty of MERCHANTABILITY or FITNESS FOR
11 A PARTICULAR PURPOSE.
12 
13 See LICENSE.txt or http://www.mitk.org for details.
14 
15 ===================================================================*/
16 
18 #include <mitkExceptionMacro.h>
19 #include <mitkSurface.h>
20 #include <vtkCell.h>
21 #include <vtkCellLinks.h>
22 #include <vtkPlane.h>
23 #include <vtkPointData.h>
24 #include <vtkPolyDataNormals.h>
25 #include <vtkSmartPointer.h>
26 
27 // forward declarations of private functions
28 static vtkIdList *GetNeighboursOfPoint(unsigned int index, vtkPolyData *polydata);
29 
30 static vtkIdList *CalculatePCAonPointNeighboursForNormalVector(int index,
31  double normal[3],
32  itk::Matrix<double, 3, 3> &mat,
33  double curVertex[3],
34  std::vector<mitk::Point3D> &pointList,
35  vtkPolyData *polyData);
36 
37 static itk::Matrix<double, 3, 3> ComputeCovarianceMatrix(itk::Matrix<double, 3, 3> &axes,
38  double sigma[3],
39  double normalizationValue);
40 namespace mitk
41 {
43  struct CovarianceMatrixCalculatorData
44  {
45  vtkPolyDataNormals *m_PolyDataNormals;
46  vtkPolyData *m_PolyData;
47  Surface *m_Input;
48  double m_VoronoiScalingFactor;
49  bool m_EnableNormalization;
50  double m_MeanVariance;
51 
52  CovarianceMatrixCalculatorData()
53  : m_PolyDataNormals(vtkPolyDataNormals::New()),
54  m_PolyData(nullptr),
55  m_Input(nullptr),
56  m_VoronoiScalingFactor(1.0),
57  m_EnableNormalization(false),
58  m_MeanVariance(0.0)
59  {
60  m_PolyDataNormals->SplittingOff();
61  }
62 
63  ~CovarianceMatrixCalculatorData()
64  {
65  if (m_PolyDataNormals)
66  m_PolyDataNormals->Delete();
67  }
68  };
69 }
70 
71 mitk::CovarianceMatrixCalculator::CovarianceMatrixCalculator() : d(new CovarianceMatrixCalculatorData())
72 {
73 }
74 
76 {
77  delete d;
78 }
79 
81 {
82  d->m_VoronoiScalingFactor = factor;
83 }
84 
86 {
87  d->m_EnableNormalization = state;
88 }
89 
91 {
92  return d->m_MeanVariance;
93 }
94 
96  const
97 {
98  return m_CovarianceMatrixList;
99 }
100 
102 {
103  d->m_Input = input;
104 }
105 
107 {
108  double normalizationValue = -1.0;
109  vtkDataArray *normals = nullptr;
110  d->m_MeanVariance = 0.0;
111 
112  if (!d->m_Input)
113  mitkThrow() << "No input surface was set in mitk::CovarianceMatrixCalculator";
114 
115  d->m_PolyData = d->m_Input->GetVtkPolyData();
116 
117  // Optional normal calculation can be disabled to use the normals
118  // of the surface:
119  // normals = d->m_PolyData->GetPointData()->GetNormals();
121  // if ( normals == NULL )
122  //{
123  d->m_PolyDataNormals->SetInputData(d->m_PolyData);
124  d->m_PolyDataNormals->Update();
125  normals = d->m_PolyDataNormals->GetOutput()->GetPointData()->GetNormals();
126  //}
127 
128  if (d->m_EnableNormalization)
129  normalizationValue = 1.5;
130 
131  // clear the matrixlist
132  m_CovarianceMatrixList.clear();
133  // allocate memory if required
134  if (d->m_PolyData->GetNumberOfPoints() > (vtkIdType)m_CovarianceMatrixList.capacity())
135  m_CovarianceMatrixList.reserve(d->m_PolyData->GetNumberOfPoints());
136 
137  for (vtkIdType i = 0; i < d->m_PolyData->GetNumberOfPoints(); ++i)
138  {
139  Vertex normal;
140  Vertex currentVertex;
141  Vertex variances = {0.0, 0.0, 0.0};
142  CovarianceMatrix mat;
143  mat.Fill(0.0);
144 
145  normals->GetTuple(i, normal);
146  d->m_PolyData->GetPoint(i, currentVertex);
147 
148  ComputeOrthonormalCoordinateSystem(i, normal, mat, variances, currentVertex);
149 
150  // use prefactor for sigma along surface
151  variances[0] = (d->m_VoronoiScalingFactor * variances[0]);
152  variances[1] = (d->m_VoronoiScalingFactor * variances[1]);
153  variances[2] = (d->m_VoronoiScalingFactor * variances[2]);
154 
155  d->m_MeanVariance += (variances[0] + variances[1] + variances[2]);
156  // compute the covariance matrix and save it
157  const CovarianceMatrix covarianceMatrix = ComputeCovarianceMatrix(mat, variances, normalizationValue);
158  m_CovarianceMatrixList.push_back(covarianceMatrix);
159  }
160 
161  if (d->m_EnableNormalization)
162  d->m_MeanVariance = normalizationValue / 3.0;
163  else
164  d->m_MeanVariance /= (3.0 * (double)d->m_PolyData->GetNumberOfPoints());
165 
166  // reset input
167  d->m_PolyData = nullptr;
168  d->m_Input = nullptr;
169 }
170 
171 // Get a list with the id's of all surrounding conected vertices
172 // to the current vertex at the given index in the polydata
173 vtkIdList *GetNeighboursOfPoint(unsigned int index, vtkPolyData *polydata)
174 {
175  vtkIdList *cellIds = vtkIdList::New();
176  vtkIdList *result = vtkIdList::New();
177  polydata->GetPointCells(index, cellIds);
178 
179  for (vtkIdType j = 0; j < cellIds->GetNumberOfIds(); j++)
180  {
181  vtkIdList *newPoints = polydata->GetCell(cellIds->GetId(j))->GetPointIds();
182  for (vtkIdType k = 0; k < newPoints->GetNumberOfIds(); k++)
183  {
184  // if point has not yet been inserted add id
185  if (result->IsId(newPoints->GetId(k)) == -1)
186  {
187  result->InsertNextId(newPoints->GetId(k));
188  }
189  }
190  }
191  cellIds->Delete();
192  return result;
193 }
194 
195 // Computes a primary component analysis of the surounding vertices
196 // of the verex at the current index.
198  double normal[3],
199  itk::Matrix<double, 3, 3> &mat,
200  double curVertex[3],
201  std::vector<mitk::Point3D> &pointList,
202  vtkPolyData *polyData)
203 {
204  typedef std::vector<mitk::Point3D> VectorType;
205  typedef VectorType::const_iterator ConstPointIterator;
206  typedef double Vertex[3];
207 
208  Vertex mean = {0.0, 0.0, 0.0};
209  Vertex tmp = {0.0, 0.0, 0.0};
210  vtkIdList *neighbourPoints = GetNeighboursOfPoint(index, polyData);
211  const vtkIdType size = neighbourPoints->GetNumberOfIds();
212 
213  // reserve memory for all neighbours
214  pointList.reserve(size);
215 
216  // project neighbours on plane given by normal
217  // and compute mean
218  for (vtkIdType i = 0; i < size; ++i)
219  {
220  mitk::Point3D p;
221  Vertex resultPoint;
222 
223  polyData->GetPoint((neighbourPoints->GetId(i)), tmp);
224 
225  vtkPlane::GeneralizedProjectPoint(tmp, curVertex, normal, resultPoint);
226 
227  p[0] = resultPoint[0];
228  p[1] = resultPoint[1];
229  p[2] = resultPoint[2];
230 
231  mean[0] += p[0];
232  mean[1] += p[1];
233  mean[2] += p[2];
234 
235  pointList.push_back(p);
236  }
237 
238  mean[0] /= (double)size;
239  mean[1] /= (double)size;
240  mean[2] /= (double)size;
241 
242  // compute the covariances with matrix multiplication
243  for (ConstPointIterator it = pointList.begin(); it != pointList.end(); ++it)
244  {
245  tmp[0] = ((*it)[0] - mean[0]);
246  tmp[1] = ((*it)[1] - mean[1]);
247  tmp[2] = ((*it)[2] - mean[2]);
248 
249  // on diagonal elements
250  mat[0][0] += tmp[0] * tmp[0];
251  mat[1][1] += tmp[1] * tmp[1];
252  mat[2][2] += tmp[2] * tmp[2];
253 
254  // of diagonal elements
255  mat[1][0] += tmp[0] * tmp[1];
256  mat[2][0] += tmp[0] * tmp[2];
257  mat[2][1] += tmp[1] * tmp[2];
258  }
259 
260  // copy upper triangle to lower triangle,
261  // we got a symetric matrix
262  mat[0][1] = mat[1][0];
263  mat[0][2] = mat[2][0];
264  mat[1][2] = mat[2][1];
265 
266  // variance
267  mat /= (size - 1);
268 
269  vnl_svd<double> svd(mat.GetVnlMatrix());
270 
271  for (int i = 0; i < 3; ++i)
272  for (int j = 0; j < 3; ++j)
273  mat[i][j] = svd.U()[j][i];
274 
275  return neighbourPoints;
276 }
277 
278 // Computes an orthonormal system for a vertex with it's surrounding neighbours.
280  const int index, Vertex normal, CovarianceMatrix &axes, Vertex variances, Vertex curVertex)
281 {
282  typedef std::vector<mitk::Point3D> VectorType;
283  typedef VectorType::const_iterator ConstPointIterator;
284  VectorType projectedPoints;
285 
286  Vertex meanValues = {0.0, 0.0, 0.0};
287  // project neighbours to new coordinate system and get principal axes
288  vtkIdList *neighbourPoints =
289  CalculatePCAonPointNeighboursForNormalVector(index, normal, axes, curVertex, projectedPoints, d->m_PolyData);
290 
291  // Set the normal as the third principal axis
292  axes[2][0] = normal[0];
293  axes[2][1] = normal[1];
294  axes[2][2] = normal[2];
295 
296  for (vtkIdType i = 0; i < neighbourPoints->GetNumberOfIds(); ++i)
297  {
298  mitk::Point3D projectedPoint;
299  Vertex curNeighbour;
300  d->m_PolyData->GetPoint(neighbourPoints->GetId(i), curNeighbour);
301 
302  curNeighbour[0] = curNeighbour[0] - curVertex[0];
303  curNeighbour[1] = curNeighbour[1] - curVertex[1];
304  curNeighbour[2] = curNeighbour[2] - curVertex[2];
305 
306  for (int k = 0; k < 3; ++k)
307  {
308  projectedPoint[k] = axes[k][0] * curNeighbour[0] + axes[k][1] * curNeighbour[1] + axes[k][2] * curNeighbour[2];
309 
310  meanValues[k] += projectedPoint[k];
311  }
312  // reuse the allocated vector from the PCA on the point neighbours
313  projectedPoints[i] = projectedPoint;
314  }
315 
316  meanValues[0] /= (double)projectedPoints.size();
317  meanValues[1] /= (double)projectedPoints.size();
318  meanValues[2] /= (double)projectedPoints.size();
319 
320  // compute variances along new axes
321  for (ConstPointIterator it = projectedPoints.begin(); it != projectedPoints.end(); ++it)
322  {
323  const mitk::Point3D &p = *it;
324 
325  variances[0] += (p[0] - meanValues[0]) * (p[0] - meanValues[0]);
326  variances[1] += (p[1] - meanValues[1]) * (p[1] - meanValues[1]);
327  variances[2] += (p[2] - meanValues[2]) * (p[2] - meanValues[2]);
328  }
329 
330  variances[0] /= (double)(projectedPoints.size() - 1);
331  variances[1] /= (double)(projectedPoints.size() - 1);
332  variances[2] /= (double)(projectedPoints.size() - 1);
333 
334  // clean up
335  neighbourPoints->Delete();
336 }
337 
338 // Sorts the axes of the computed orthonormal system based on
339 // the eigenvalues in a descending order
340 itk::Matrix<double, 3, 3> ComputeCovarianceMatrix(itk::Matrix<double, 3, 3> &axes,
341  double sigma[3],
342  double normalizationValue)
343 {
344  unsigned int idxMax, idxMin, idxBetween;
345  itk::Matrix<double, 3, 3> returnValue;
346  itk::Matrix<double, 3, 3> V;
347  itk::Matrix<double, 3, 3> diagMatrix;
348  diagMatrix.Fill(0.0);
349 
350  if (sigma[0] >= sigma[1] && sigma[0] >= sigma[2])
351  {
352  idxMax = 0;
353  if (sigma[1] >= sigma[2])
354  {
355  idxBetween = 1;
356  idxMin = 2;
357  }
358  else
359  {
360  idxBetween = 2;
361  idxMin = 1;
362  }
363  }
364  else if (sigma[1] >= sigma[0] && sigma[1] >= sigma[2])
365  {
366  idxMax = 1;
367  if (sigma[0] >= sigma[2])
368  {
369  idxBetween = 0;
370  idxMin = 2;
371  }
372  else
373  {
374  idxBetween = 2;
375  idxMin = 0;
376  }
377  }
378  else // index 2 corresponds to largest sigma
379  {
380  idxMax = 2;
381  if (sigma[0] >= sigma[1])
382  {
383  idxBetween = 0;
384  idxMin = 1;
385  }
386  else
387  {
388  idxBetween = 1;
389  idxMin = 0;
390  }
391  }
392 
393  V[0][0] = axes[idxMax][0];
394  V[1][0] = axes[idxMax][1];
395  V[2][0] = axes[idxMax][2];
396  V[0][1] = axes[idxBetween][0];
397  V[1][1] = axes[idxBetween][1];
398  V[2][1] = axes[idxBetween][2];
399  V[0][2] = axes[idxMin][0];
400  V[1][2] = axes[idxMin][1];
401  V[2][2] = axes[idxMin][2];
402 
403  diagMatrix[0][0] = sigma[idxMax];
404  diagMatrix[1][1] = sigma[idxBetween];
405  diagMatrix[2][2] = sigma[idxMin];
406 
407  returnValue = V * diagMatrix * V.GetTranspose();
408 
409  if (normalizationValue > 0.0)
410  {
411  double trace = returnValue[0][0] + returnValue[1][1] + returnValue[2][2];
412  returnValue *= (normalizationValue / trace);
413  }
414 
415  return returnValue;
416 }
Class for storing surfaces (vtkPolyData).
Definition: mitkSurface.h:32
std::vector< CovarianceMatrix > CovarianceMatrixList
DataCollection - Class to facilitate loading/accessing structured data.
const CovarianceMatrixList & GetCovarianceMatrices() const
itk::Vector< float, 3 > VectorType
static vtkIdList * CalculatePCAonPointNeighboursForNormalVector(int index, double normal[3], itk::Matrix< double, 3, 3 > &mat, double curVertex[3], std::vector< mitk::Point3D > &pointList, vtkPolyData *polyData)
static vtkIdList * GetNeighboursOfPoint(unsigned int index, vtkPolyData *polydata)
void ComputeOrthonormalCoordinateSystem(const int index, Vertex normal, CovarianceMatrix &principalComponents, Vertex variances, Vertex curVertex)
#define mitkThrow()
static itk::Matrix< double, 3, 3 > ComputeCovarianceMatrix(itk::Matrix< double, 3, 3 > &axes, double sigma[3], double normalizationValue)
static itkEventMacro(BoundingShapeInteractionEvent, itk::AnyEvent) class MITKBOUNDINGSHAPE_EXPORT BoundingShapeInteractor Pointer New()
Basic interaction methods for mitk::GeometryData.