Chaste  Build::
VesselNetworkReader.cpp
1 /*
2 
3  Copyright (c) 2005-2015, University of Oxford.
4  All rights reserved.
5 
6  University of Oxford means the Chancellor, Masters and Scholars of the
7  University of Oxford, having an administrative office at Wellington
8  Square, Oxford OX1 2JD, UK.
9 
10  This file is part of Chaste.
11 
12  Redistribution and use in source and binary forms, with or without
13  modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24  AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25  IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26  ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27  LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28  CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29  GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30  HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31  LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32  OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34  */
35 
36 #define _BACKWARD_BACKWARD_WARNING_H 1 //Cut out the strstream deprecated warning for now (gcc4.3)
37 #include <vtkDoubleArray.h>
38 #include <vtkCellData.h>
39 #include <vtkCellArray.h>
40 #include <vtkPointData.h>
41 #include <vtkPolyData.h>
42 #include <vtkLine.h>
43 #include <vtkXMLPolyDataReader.h>
44 #include <vtkSmartPointer.h>
45 #include <vtkCleanPolyData.h>
46 #include <vtkSplineFilter.h>
47 #include <vtkVersion.h>
48 #include "Exception.hpp"
49 #include "VesselNetworkReader.hpp"
50 
51 template<unsigned DIM>
53  : mFileName(),
54  mRadiusLabel("Node Radius"),
55  mRadiusConversionFactor(1.0),
56  mMergeCoincidentPoints(false),
57  mTargetSegmentLength(0.0*unit::metres),
58  mReferenceLength(1.e-6 * unit::metres)
59 {
60 
61 }
62 
63 template<unsigned DIM>
65 {
66 }
67 
68 template <unsigned DIM>
69 boost::shared_ptr<VesselNetworkReader<DIM> > VesselNetworkReader<DIM>::Create()
70 {
71  MAKE_PTR(VesselNetworkReader<DIM>, pSelf);
72  return pSelf;
73 }
74 
75 template <unsigned DIM>
76 void VesselNetworkReader<DIM>::SetRadiusArrayName(const std::string& rRadius)
77 {
78  mRadiusLabel = rRadius;
79 }
80 
81 template <unsigned DIM>
83 {
84  mMergeCoincidentPoints = mergePoints;
85 }
86 
87 template <unsigned DIM>
88 void VesselNetworkReader<DIM>::SetTargetSegmentLength(units::quantity<unit::length> targetSegmentLength)
89 {
90  mTargetSegmentLength = targetSegmentLength;
91 }
92 
93 template <unsigned DIM>
94 void VesselNetworkReader<DIM>::SetReferenceLengthScale(units::quantity<unit::length> rReferenceLength)
95 {
96  mReferenceLength = rReferenceLength;
97 }
98 
99 template<unsigned DIM>
100 boost::shared_ptr<VesselNetwork<DIM> > VesselNetworkReader<DIM>::Read()
101 {
102  if(mFileName.empty())
103  {
104  EXCEPTION("File name not set in vessel network reader");
105  }
106 
107  // Create an empty vessel network
108  boost::shared_ptr<VesselNetwork<DIM> > p_network = VesselNetwork<DIM>::Create();
109 
110  // Create a VTK PolyData object based on the contents of the input VTK file
111  vtkSmartPointer<vtkXMLPolyDataReader> p_reader = vtkSmartPointer<vtkXMLPolyDataReader>::New();
112  p_reader->SetFileName(mFileName.c_str());
113  p_reader->Update();
114 
115  vtkSmartPointer<vtkPolyData> p_polydata = p_reader->GetOutput();
116 
118  {
119  vtkSmartPointer<vtkCleanPolyData> p_cleaner = vtkSmartPointer<vtkCleanPolyData>::New();
120  #if VTK_MAJOR_VERSION <= 5
121  p_cleaner->SetInput(p_polydata);
122  #else
123  p_cleaner->SetInputData(p_polydata);
124  #endif
125  p_cleaner->Update();
126  p_polydata = p_cleaner->GetOutput();
127  }
128 
129  if(mTargetSegmentLength != 0.0*unit::metres)
130  {
131  vtkSmartPointer<vtkSplineFilter> p_spline_filter = vtkSmartPointer<vtkSplineFilter>::New();
132  #if VTK_MAJOR_VERSION <= 5
133  p_spline_filter->SetInput(p_polydata);
134  #else
135  p_spline_filter->SetInputData(p_polydata);
136  #endif
137  p_spline_filter->SetSubdivideToLength();
138  p_spline_filter->SetLength(mTargetSegmentLength/mReferenceLength);
139  p_spline_filter->Update();
140  p_polydata = p_spline_filter->GetOutput();
141  }
142 
143  // Create the nodes
144  vtkSmartPointer<vtkPointData> p_point_data = vtkSmartPointer<vtkPointData>::New();
145  std::vector<boost::shared_ptr<VesselNode<DIM> > > nodes;
146  for (vtkIdType i = 0; i < p_polydata->GetNumberOfPoints(); i++)
147  {
148  double point_coords[3];
149  p_polydata->GetPoint(i, point_coords);
150  if (DIM < 3)
151  {
152  nodes.push_back(VesselNode<DIM>::Create(point_coords[0], point_coords[1], 0.0, mReferenceLength));
153  }
154  else
155  {
156  nodes.push_back(VesselNode<DIM>::Create(point_coords[0], point_coords[1], point_coords[2], mReferenceLength));
157  }
158  }
159 
160  // Extract radii corresponding to each node from the VTK Polydata and store them in a list.
161  p_point_data = p_polydata->GetPointData();
162  std::vector<double> radii;
163  for (vtkIdType i = 0; i < p_point_data->GetNumberOfArrays(); i++)
164  {
165  std::string array_name = p_point_data->GetArrayName(i);
166  if (array_name.compare(mRadiusLabel) == 0)
167  {
168  for (vtkIdType j = 0; j < p_point_data->GetArray(i)->GetNumberOfTuples(); j++)
169  {
170  radii.push_back(p_point_data->GetArray(i)->GetTuple1(j));
171  }
172  }
173  }
174 
175  // Extract vessels from the VTK Polydata. This is done by iterating over a VTK CellArray object which
176  // returns a 'pointList' vtkIdList object. This object contains the point IDs of the nodes which make up
177  // the vessel.
178  vtkSmartPointer<vtkCellArray> pCellArray = vtkSmartPointer<vtkCellArray>::New();
179  pCellArray = p_polydata->GetLines();
180  for (int i = 0; i < p_polydata->GetNumberOfLines(); i++)
181  {
182  // Make a new vessel
183  vtkIdType num_segments;
184  vtkIdType* pSegmentList;
185  pCellArray->GetNextCell(num_segments, pSegmentList);
186  std::vector<boost::shared_ptr<VesselSegment<DIM> > > segments;
187  // Add segments to the vessels in order
188  for (int j = 1; j < num_segments; j++)
189  {
190  boost::shared_ptr<VesselSegment<DIM> > p_segment = VesselSegment<DIM>::Create(nodes[pSegmentList[j - 1]],nodes[pSegmentList[j]]);
191  if(unsigned(radii.size())> pSegmentList[j])
192  {
193  p_segment->SetRadius(radii[pSegmentList[j]] * mReferenceLength);
194  }
195  segments.push_back(p_segment);
196  }
197  // Add the resulting vessel to the network
198  p_network->AddVessel(Vessel<DIM>::Create(segments));
199  }
200  return p_network;
201 }
202 
203 template<unsigned DIM>
204 void VesselNetworkReader<DIM>::SetFileName(const std::string& rFileName)
205 {
206  mFileName = rFileName;
207 }
208 
209 //Explicit instantiation
210 template class VesselNetworkReader<2> ;
211 template class VesselNetworkReader<3> ;
void SetMergeCoincidentPoints(bool mergePoints)
static boost::shared_ptr< VesselNetworkReader< DIM > > Create()
void SetReferenceLengthScale(units::quantity< unit::length > rReferenceLength)
units::quantity< unit::length > mReferenceLength
static boost::shared_ptr< VesselNetwork< DIM > > Create()
units::quantity< unit::length > mTargetSegmentLength
boost::shared_ptr< VesselNetwork< DIM > > Read()
static boost::shared_ptr< VesselSegment< DIM > > Create(boost::shared_ptr< VesselNode< DIM > > pNode1, boost::shared_ptr< VesselNode< DIM > > pNode2)
void SetRadiusArrayName(const std::string &rRadius)
void SetTargetSegmentLength(units::quantity< unit::length > targetSegmentLength)
void SetFileName(const std::string &rFileName)