CustusX  18.04
An IGT application
cxBranchList.cpp
Go to the documentation of this file.
1 /*=========================================================================
2 This file is part of CustusX, an Image Guided Therapy Application.
3 
4 Copyright (c) SINTEF Department of Medical Technology.
5 All rights reserved.
6 
7 CustusX is released under a BSD 3-Clause license.
8 
9 See Lisence.txt (https://github.com/SINTEFMedtek/CustusX/blob/master/License.txt) for details.
10 =========================================================================*/
11 #include "cxBranchList.h"
12 #include "cxBranch.h"
13 #include "cxMesh.h"
14 #include "cxVector3D.h"
15 #include <vtkPolyData.h>
16 #include <vtkCardinalSpline.h>
17 
18 
19 typedef vtkSmartPointer<class vtkCardinalSpline> vtkCardinalSplinePtr;
20 
21 namespace cx
22 {
23 
25 {
26 
27 }
28 
29 
31 {
32 // for (int i = 0; i < mBranches.size(); i++)
33 // mBranches[i]->~Branch();
34 }
35 
37 {
38  mBranches.push_back(b);
39 }
40 
42 {
43  for( int i = 0; i < mBranches.size(); i++ )
44  {
45  if (b == mBranches[i])
46  {
47  mBranches.erase(mBranches.begin() + i);
48  return;
49  }
50  }
51 }
52 
54 {
55  mBranches.clear();
56 }
57 
58 std::vector<BranchPtr> BranchList::getBranches()
59 {
60  return mBranches;
61 }
62 
63 void BranchList::selectGenerations(int maxGeneration)
64 {
65  std::vector<int> branchNumbersToBeDeleted;
66  for( int i = 0; i < mBranches.size(); i++ )
67  {
68  int generationCounter = 1;
69  BranchPtr currentBranch = mBranches[i];
70  while (currentBranch->getParentBranch()){
71  generationCounter++;
72  currentBranch = currentBranch->getParentBranch();
73  if (generationCounter > maxGeneration)
74  {
75  branchNumbersToBeDeleted.push_back(i);
76  break;
77  }
78  }
79 
80  }
81 
82  for ( int i = branchNumbersToBeDeleted.size() - 1; i >= 0; i-- )
83  deleteBranch(mBranches[branchNumbersToBeDeleted[i]]);
84 }
85 
86 
88 {
89  for (int i = 0; i < mBranches.size(); i++)
90  {
91  Eigen::MatrixXd positions = mBranches[i]->getPositions();
92  Eigen::MatrixXd diff = positions.rightCols(positions.cols() - 1) - positions.leftCols(positions.cols() - 1);
93  Eigen::MatrixXd orientations(positions.rows(),positions.cols());
94  orientations.leftCols(orientations.cols() - 1) = diff / diff.norm();
95  orientations.rightCols(1) = orientations.col(orientations.cols() - 1);
96  mBranches[i]->setOrientations(orientations);
97  }
98 }
99 
101 {
102  for (int i = 0; i < mBranches.size(); i++)
103  {
104  Eigen::MatrixXd orientations = mBranches[i]->getOrientations();
105  Eigen::MatrixXd newOrientations(orientations.rows(),orientations.cols());
106  int numberOfColumns = orientations.cols();
107  for (int j = 0; j < numberOfColumns; j++)
108  {
109  newOrientations.col(j) = orientations.block(0,std::max(j-2,0),orientations.rows(),std::min(5,numberOfColumns-j)).rowwise().mean(); //smoothing
110  newOrientations.col(j) = newOrientations.col(j) / newOrientations.col(j).norm(); // normalizing
111  }
112  mBranches[i]->setOrientations(newOrientations);
113  }
114 }
115 
116 void BranchList::interpolateBranchPositions(int interpolationFactor){
117 
118  for (int i = 0; i < mBranches.size(); i++)
119  {
120  Eigen::MatrixXd positions = mBranches[i]->getPositions();
121  std::vector<Eigen::Vector3d> interpolatedPositions;
122  for (int j = 0; j < positions.cols()-1; j++)
123  {
124  for (int k = 0; k < interpolationFactor; k++){
125  Eigen::Vector3d interpolationPoint;
126  interpolationPoint[0] = (positions(0,j)*(interpolationFactor-k) + positions(0,j+1)*(k) ) / interpolationFactor;
127  interpolationPoint[1] = (positions(1,j)*(interpolationFactor-k) + positions(1,j+1)*(k) ) / interpolationFactor;
128  interpolationPoint[2] = (positions(2,j)*(interpolationFactor-k) + positions(2,j+1)*(k) ) / interpolationFactor;
129  interpolatedPositions.push_back(interpolationPoint);
130  }
131  }
132  Eigen::MatrixXd interpolationResult(3 , interpolatedPositions.size());
133  for (int j = 0; j < interpolatedPositions.size(); j++)
134  {
135  interpolationResult(0,j) = interpolatedPositions[j](0);
136  interpolationResult(1,j) = interpolatedPositions[j](1);
137  interpolationResult(2,j) = interpolatedPositions[j](2);
138  }
139  mBranches[i]->setPositions(interpolationResult);
140  }
141 
142 }
143 
144 void BranchList::smoothBranchPositions(int controlPointDistance)
145 {
146  for (int i = 0; i < mBranches.size(); i++)
147  {
148  Eigen::MatrixXd positions = mBranches[i]->getPositions();
149  int numberOfInputPoints = positions.cols();
150  //int controlPointFactor = 10;
151  //int numberOfControlPoints = numberOfInputPoints / controlPointFactor;
152  double branchLength = (positions.rightCols(1) - positions.leftCols(1)).norm();
153  int numberOfControlPoints = std::ceil(branchLength/controlPointDistance);
154  numberOfControlPoints = std::max(numberOfControlPoints, 2); // at least two control points
155 
156  vtkCardinalSplinePtr splineX = vtkSmartPointer<vtkCardinalSpline>::New();
157  vtkCardinalSplinePtr splineY = vtkSmartPointer<vtkCardinalSpline>::New();
158  vtkCardinalSplinePtr splineZ = vtkSmartPointer<vtkCardinalSpline>::New();
159 
160  //add control points to spline
161  for(int j=0; j<numberOfControlPoints; j++)
162  {
163  int indexP = (j*numberOfInputPoints)/numberOfControlPoints;
164 
165  splineX->AddPoint(indexP,positions(0,indexP));
166  splineY->AddPoint(indexP,positions(1,indexP));
167  splineZ->AddPoint(indexP,positions(2,indexP));
168  }
169  //Always add the last point to complete spline
170  splineX->AddPoint(numberOfInputPoints-1,positions(0,numberOfInputPoints-1));
171  splineY->AddPoint(numberOfInputPoints-1,positions(1,numberOfInputPoints-1));
172  splineZ->AddPoint(numberOfInputPoints-1,positions(2,numberOfInputPoints-1));
173 
174  //evaluate spline - get smoothed positions
175  Eigen::MatrixXd smoothingResult(3 , numberOfInputPoints);
176  for(int j=0; j<numberOfInputPoints; j++)
177  {
178  double splineParameter = j;
179  smoothingResult(0,j) = splineX->Evaluate(splineParameter);
180  smoothingResult(1,j) = splineY->Evaluate(splineParameter);
181  smoothingResult(2,j) = splineZ->Evaluate(splineParameter);
182  }
183  mBranches[i]->setPositions(smoothingResult);
184  }
185 }
186 
187 void BranchList::findBranchesInCenterline(Eigen::MatrixXd positions_r)
188 {
189  positions_r = sortMatrix(2,positions_r);
190  Eigen::MatrixXd positionsNotUsed_r = positions_r;
191 
192 // int minIndex;
193  int index;
194  int splitIndex;
195  Eigen::MatrixXd::Index startIndex;
196  BranchPtr branchToSplit;
197  while (positionsNotUsed_r.cols() > 0)
198  {
199  if (!mBranches.empty())
200  {
201  double minDistance = 1000;
202  for (int i = 0; i < mBranches.size(); i++)
203  {
204  std::pair<std::vector<Eigen::MatrixXd::Index>, Eigen::VectorXd> distances;
205  distances = dsearchn(positionsNotUsed_r, mBranches[i]->getPositions());
206  double d = distances.second.minCoeff(&index);
207  if (d < minDistance)
208  {
209  minDistance = d;
210  branchToSplit = mBranches[i];
211  startIndex = index;
212  if (minDistance < 2)
213  break;
214  }
215  }
216  std::pair<Eigen::MatrixXd::Index, double> dsearchResult = dsearch(positionsNotUsed_r.col(startIndex) , branchToSplit->getPositions());
217  splitIndex = dsearchResult.first;
218  }
219  else //if this is the first branch. Select the top position (Trachea).
220  startIndex = positionsNotUsed_r.cols() - 1;
221 
222  std::pair<Eigen::MatrixXd,Eigen::MatrixXd > connectedPointsResult = findConnectedPointsInCT(startIndex , positionsNotUsed_r);
223  Eigen::MatrixXd branchPositions = connectedPointsResult.first;
224  positionsNotUsed_r = connectedPointsResult.second;
225 
226  if (branchPositions.cols() >= 5) //only include brances of length >= 5 points
227  {
228  BranchPtr newBranch = BranchPtr(new Branch());
229  newBranch->setPositions(branchPositions);
230  mBranches.push_back(newBranch);
231 
232  if (mBranches.size() > 1) // do not try to split another branch when the first branch is processed
233  {
234  if ((splitIndex + 1 >= 5) && (branchToSplit->getPositions().cols() - splitIndex - 1 >= 5))
235  //do not split branch if the new branch is close to the edge of the branch
236  //if the new branch is not close to one of the edges of the
237  //connected existing branch: Split the existing branch
238  {
239  BranchPtr newBranchFromSplit = BranchPtr(new Branch());
240  Eigen::MatrixXd branchToSplitPositions = branchToSplit->getPositions();
241  newBranchFromSplit->setPositions(branchToSplitPositions.rightCols(branchToSplitPositions.cols() - splitIndex - 1));
242  branchToSplit->setPositions(branchToSplitPositions.leftCols(splitIndex + 1));
243  mBranches.push_back(newBranchFromSplit);
244  newBranchFromSplit->setParentBranch(branchToSplit);
245  newBranch->setParentBranch(branchToSplit);
246  newBranchFromSplit->setChildBranches(branchToSplit->getChildBranches());
247  branchVector branchToSplitChildren = branchToSplit->getChildBranches();
248  for (int i = 0; i < branchToSplitChildren.size(); i++)
249  branchToSplitChildren[i]->setParentBranch(newBranchFromSplit);
250  branchToSplit->deleteChildBranches();
251  branchToSplit->addChildBranch(newBranchFromSplit);
252  branchToSplit->addChildBranch(newBranch);
253  }
254  else if (splitIndex + 1 < 5)
255  // If the new branch is close to the start of the existing
256  // branch: Connect it to the same position start as the
257  // existing branch
258  {
259  newBranch->setParentBranch(branchToSplit->getParentBranch());
260  if(branchToSplit->getParentBranch())
261  branchToSplit->getParentBranch()->addChildBranch(newBranch);
262  }
263  else if (branchToSplit->getPositions().cols() - splitIndex - 1 < 5)
264  // If the new branch is close to the end of the existing
265  // branch: Connect it to the end of the existing branch
266  {
267  newBranch->setParentBranch(branchToSplit);
268  branchToSplit->addChildBranch(newBranch);
269  }
270 
271  }
272 
273  }
274  }
275 }
276 
277 
278 BranchListPtr BranchList::removePositionsForLocalRegistration(Eigen::MatrixXd trackingPositions, double maxDistance)
279 {
280  BranchListPtr retval = BranchListPtr(new BranchList());
281  BranchPtr b;
282  for (int i = 0; i < mBranches.size(); i++)
283  {
284  b = BranchPtr(new Branch());
285  b->setPositions(mBranches[i]->getPositions());
286  b->setOrientations(mBranches[i]->getOrientations());
287  retval->addBranch(b);
288  }
289 
290  std::vector<BranchPtr> branches = retval->getBranches();
291  Eigen::MatrixXd positions;
292  Eigen::MatrixXd orientations;
293  for (int i = 0; i < branches.size(); i++)
294  {
295  positions = branches[i]->getPositions();
296  orientations = branches[i]->getOrientations();
297  std::pair<std::vector<Eigen::MatrixXd::Index>, Eigen::VectorXd> distanceData;
298  distanceData = dsearchn(positions, trackingPositions);
299  Eigen::VectorXd distance = distanceData.second;
300  for (int j = positions.cols() - 1; j >= 0; j--)
301  {
302  if (distance(j) > maxDistance)
303  {
304  positions = eraseCol(j, positions);
305  orientations = eraseCol(j, orientations);
306  }
307  }
308  branches[i]->setPositions(positions);
309  branches[i]->setOrientations(orientations);
310  }
311  return retval;
312 }
313 
332 vtkPolyDataPtr BranchList::createVtkPolyDataFromBranches(bool fullyConnected, bool straightBranches) const
333 {
334  vtkPolyDataPtr retval = vtkPolyDataPtr::New();
335  vtkPointsPtr points = vtkPointsPtr::New();
336  vtkCellArrayPtr lines = vtkCellArrayPtr::New();
337 
338  int positionCounter = 0;
339  for (size_t i = 0; i < mBranches.size(); ++i)
340  {
341  Eigen::MatrixXd positions = mBranches[i]->getPositions();
342  if(straightBranches)
343  {
344  ++positionCounter;
345  points->InsertNextPoint(positions(0,0),positions(1,0),positions(2,0));
346  points->InsertNextPoint(positions(0,positions.cols()-1),positions(1,positions.cols()-1),positions(2,positions.cols()-1));
347  vtkIdType connection[2] = {positionCounter-1, positionCounter};
348  lines->InsertNextCell(2, connection);
349  ++positionCounter;
350  }
351  else
352  {
353  for (int j = 0; j < positions.cols(); ++j)
354  {
355  ++positionCounter;
356  points->InsertNextPoint(positions(0,j),positions(1,j),positions(2,j));
357  if (j < positions.cols()-1)
358  {
359  vtkIdType connection[2] = {positionCounter-1, positionCounter};
360  lines->InsertNextCell(2, connection);
361  }
362  }
363  }
364  }
365  if(fullyConnected)
366  {
367  int this_branchs_first_point_in_full_polydata_point_list = 0;
368  for(size_t i = 0; i < mBranches.size(); ++i)
369  {
370  if(i>0)
371  {
372  if(!straightBranches)
373  this_branchs_first_point_in_full_polydata_point_list += mBranches[i-1]->getPositions().cols();
374  else
375  this_branchs_first_point_in_full_polydata_point_list += 2;
376  }
377  int parent_index_in_branch_list = mBranches[i]->findParentIndex(mBranches);
378 
379  if(parent_index_in_branch_list > -1)
380  {
381  int parent_branch_last_point_in_full_polydata = 0;
382  for(int j = 0; j <= parent_index_in_branch_list; ++j)
383  {
384  if(!straightBranches)
385  parent_branch_last_point_in_full_polydata += mBranches[j]->getPositions().cols() - 1;
386  else
387  parent_branch_last_point_in_full_polydata += (1 + j*2);
388  }
389  vtkIdType connection[2] = {parent_branch_last_point_in_full_polydata, this_branchs_first_point_in_full_polydata_point_list};
390  lines->InsertNextCell(2, connection);
391  }
392 
393  }
394 
395  }
396  retval->SetPoints(points);
397  retval->SetLines(lines);
398 
399  return retval;
400 }
401 
402 Eigen::MatrixXd sortMatrix(int rowNumber, Eigen::MatrixXd matrix)
403 {
404  for (int i = 0; i < matrix.cols() - 1; i++) {
405  for (int j = i + 1; j < matrix.cols(); j++) {
406  if (matrix(rowNumber,i) > matrix(rowNumber,j)){
407  matrix.col(i).swap(matrix.col(j));
408  }
409  }
410  }
411 return matrix;
412 }
413 
414 
415 
416 Eigen::MatrixXd eraseCol(int removeIndex, Eigen::MatrixXd positions)
417 {
418  positions.block(0 , removeIndex , positions.rows() , positions.cols() - removeIndex - 1) = positions.rightCols(positions.cols() - removeIndex - 1);
419  positions.conservativeResize(Eigen::NoChange, positions.cols() - 1);
420  return positions;
421 }
422 
423 std::pair<Eigen::MatrixXd::Index, double> dsearch(Eigen::Vector3d p, Eigen::MatrixXd positions)
424 {
425  Eigen::MatrixXd::Index index;
426  // find nearest neighbour
427  (positions.colwise() - p).colwise().squaredNorm().minCoeff(&index);
428  double d = (positions.col(index) - p).norm();
429 
430  return std::make_pair(index , d);
431 }
432 
433 std::pair<std::vector<Eigen::MatrixXd::Index>, Eigen::VectorXd > dsearchn(Eigen::MatrixXd p1, Eigen::MatrixXd p2)
434 {
435  Eigen::MatrixXd::Index index;
436  std::vector<Eigen::MatrixXd::Index> indexVector;
437  Eigen::VectorXd D(p1.cols());
438  for (int i = 0; i < p1.cols(); i++)
439  {
440  // find nearest neighbour
441  (p2.colwise() - p1.col(i)).colwise().squaredNorm().minCoeff(&index);
442  D(i) = (p2.col(index) - p1.col(i)).norm();
443  indexVector.push_back(index);
444  }
445  return std::make_pair(indexVector , D);
446 }
447 
448 std::pair<Eigen::MatrixXd,Eigen::MatrixXd > findConnectedPointsInCT(int startIndex , Eigen::MatrixXd positionsNotUsed)
449 {
450  //Eigen::MatrixXd branchPositions(positionsNotUsed.rows(), positionsNotUsed.cols());
451  Eigen::MatrixXd thisPosition(3,1);
452  std::vector<Eigen::MatrixXd> branchPositionsVector;
453  thisPosition = positionsNotUsed.col(startIndex);
454  branchPositionsVector.push_back(thisPosition); //add first position to branch
455  positionsNotUsed = eraseCol(startIndex,positionsNotUsed);; //remove first position from list of remaining points
456 
457  while (positionsNotUsed.cols() > 0)
458  {
459  std::pair<Eigen::MatrixXd::Index, double > minDistance = dsearch(thisPosition, positionsNotUsed);
460  Eigen::MatrixXd::Index index = minDistance.first;
461  double d = minDistance.second;
462  if (d > 3) // more than 3 mm distance to closest point --> branch is compledted
463  break;
464 
465  thisPosition = positionsNotUsed.col(index);
466  positionsNotUsed = eraseCol(index,positionsNotUsed);
467  //add position to branch
468  branchPositionsVector.push_back(thisPosition);
469 
470  }
471 
472  Eigen::MatrixXd branchPositions(3,branchPositionsVector.size());
473 
474  for (int j = 0; j < branchPositionsVector.size(); j++)
475  {
476  branchPositions.block(0,j,3,1) = branchPositionsVector[j];
477  }
478 
479  return std::make_pair(branchPositions, positionsNotUsed);
480 }
481 
482 
483 }//namespace cx
virtual ~BranchList()
Vector3D ceil(const Vector3D &a)
Definition: cxVector3D.cpp:84
std::pair< Eigen::MatrixXd, Eigen::MatrixXd > findConnectedPointsInCT(int startIndex, Eigen::MatrixXd positionsNotUsed)
void interpolateBranchPositions(int interpolationFactor)
vtkSmartPointer< class vtkCellArray > vtkCellArrayPtr
boost::shared_ptr< class BranchList > BranchListPtr
BranchListPtr removePositionsForLocalRegistration(Eigen::MatrixXd trackingPositions, double maxDistance)
vtkSmartPointer< vtkPoints > vtkPointsPtr
void addBranch(BranchPtr b)
void deleteAllBranches()
boost::shared_ptr< class Branch > BranchPtr
vtkSmartPointer< class vtkCardinalSpline > vtkCardinalSplinePtr
vtkSmartPointer< class vtkCardinalSpline > vtkCardinalSplinePtr
Definition: cxAccusurf.cpp:17
vtkPolyDataPtr createVtkPolyDataFromBranches(bool fullyConnected=false, bool straightBranches=false) const
BranchList::createVtkPolyDataFromBranches Return a VtkPolyData object created from the branches in th...
std::pair< Eigen::MatrixXd::Index, double > dsearch(Eigen::Vector3d p, Eigen::MatrixXd positions)
void calculateOrientations()
void smoothBranchPositions(int controlPointDistance)
std::vector< BranchPtr > getBranches()
void deleteBranch(BranchPtr b)
vtkSmartPointer< vtkPolyData > vtkPolyDataPtr
void findBranchesInCenterline(Eigen::MatrixXd positions_r)
class org_custusx_registration_method_bronchoscopy_EXPORT Branch
Definition: cxBranch.h:26
std::pair< std::vector< Eigen::MatrixXd::Index >, Eigen::VectorXd > dsearchn(Eigen::MatrixXd p1, Eigen::MatrixXd p2)
Eigen::MatrixXd eraseCol(int removeIndex, Eigen::MatrixXd positions)
void selectGenerations(int maxGeneration)
Eigen::MatrixXd sortMatrix(int rowNumber, Eigen::MatrixXd matrix)
void smoothOrientations()
std::vector< BranchPtr > branchVector
Definition: cxBranch.h:28
Namespace for all CustusX production code.