Chaste  Build::
NetworkToSurface.cpp
1 /*
2 
3 Copyright (c) 2005-2016, 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 #include "NetworkToSurface.hpp"
37 #include "NetworkToImage.hpp"
38 #include <vtkMarchingSquares.h>
39 #include <vtkMarchingCubes.h>
40 #include <vtkCleanPolyData.h>
41 #include <vtkStripper.h>
42 #include <vtkSplineFilter.h>
43 #include <vtkTriangleFilter.h>
44 #include <vtkTransform.h>
45 #include <vtkClipPolyData.h>
46 #include <vtkPolyData.h>
47 #include <vtkIdList.h>
48 #include <vtkLine.h>
49 #include <vtkCellArray.h>
50 #include <vtkBox.h>
51 #include <vtkVersion.h>
52 #include <vtkProbeFilter.h>
53 #include <vtkPolyDataNormals.h>
54 #include <vtkWindowedSincPolyDataFilter.h>
55 //#include <vtkvmtkSimpleCapPolyData.h>
56 //#include <vtkvmtkPolyDataSizingFunction.h>
57 //#include <vtkvmtkPolyDataToUnstructuredGridFilter.h>
58 #include "UblasCustomFunctions.hpp"
59 #include "UblasIncludes.hpp"
60 #include "UnitCollection.hpp"
61 
62 template<unsigned DIM>
64  mpNetwork(),
65  mSamplingGridSpacing(2.0 * 1.e-6 * unit::metres),
66  mSplineResamplingLength(10.0 * 1.e-6 * unit::metres),
67  mpSurface(),
68  mpImage(),
69  mReferenceLength(BaseUnits::Instance()->GetReferenceLengthScale())
70 {
71 
72 }
73 
74 template<unsigned DIM>
75 boost::shared_ptr<NetworkToSurface<DIM> > NetworkToSurface<DIM>::Create()
76 {
77  MAKE_PTR(NetworkToSurface<DIM>, pSelf);
78  return pSelf;
79 }
80 
81 template<unsigned DIM>
83 {
84 
85 }
86 
87 template<unsigned DIM>
88 vtkSmartPointer<vtkPolyData> NetworkToSurface<DIM>::GetSurface()
89 {
90  return mpSurface;
91 }
92 
93 template<unsigned DIM>
94 vtkSmartPointer<vtkImageData> NetworkToSurface<DIM>::GetSamplingImage()
95 {
96  return mpImage;
97 }
98 
99 template<unsigned DIM>
101 {
102  mpNetwork = pNetwork;
103 }
104 
105 template<unsigned DIM>
106 void NetworkToSurface<DIM>::SetResamplingSplineSize(units::quantity<unit::length> splineResampleSize)
107 {
108  mSplineResamplingLength = splineResampleSize;
109 }
110 
111 template<unsigned DIM>
112 void NetworkToSurface<DIM>::SetResamplingGridSize(units::quantity<unit::length> sampleGridSize)
113 {
114  mSamplingGridSpacing = sampleGridSize;
115 }
116 
117 template<unsigned DIM>
119 {
120  if(DIM==2)
121  {
122  // Interpolate the network onto a regular grid
123  boost::shared_ptr<NetworkToImage<DIM> > p_net_to_image = NetworkToImage<DIM>::Create();
124  p_net_to_image->SetNetwork(mpNetwork);
125  p_net_to_image->SetGridSpacing(mSamplingGridSpacing);
126  p_net_to_image->SetPaddingFactors(0.1, 0.1, 0.0);
127  p_net_to_image->SetImageDimension(2);
128  p_net_to_image->Update();
129  mpImage = p_net_to_image->GetOutput();
130 
131  // Get the outer boundaries of the network
132  vtkSmartPointer<vtkMarchingSquares> p_squares = vtkSmartPointer<vtkMarchingSquares>::New();
133  #if VTK_MAJOR_VERSION <= 5
134  p_squares->SetInput(mpImage);
135  #else
136  p_squares->SetInputData(mpImage);
137  #endif
138  p_squares->SetValue(0, 1);
139 
140  // Remove duplicate points and join up lines to form polylines
141  vtkSmartPointer<vtkCleanPolyData> p_cleaner = vtkSmartPointer<vtkCleanPolyData>::New();
142  p_cleaner->SetInputConnection(p_squares->GetOutputPort());
143 
144  vtkSmartPointer<vtkStripper> p_stripper = vtkSmartPointer<vtkStripper>::New();
145  p_stripper->SetInputConnection(p_cleaner->GetOutputPort());
146 
147  // Downsample and smooth the polyline
148  vtkSmartPointer<vtkSplineFilter> p_spline= vtkSmartPointer<vtkSplineFilter>::New();
149  p_spline->SetLength(mSplineResamplingLength/mReferenceLength);
150  p_spline->SetSubdivideToLength();
151  p_spline->SetInputConnection(p_stripper->GetOutputPort());
152 
153  vtkSmartPointer<vtkTriangleFilter> p_triangle = vtkSmartPointer<vtkTriangleFilter>::New();
154  p_triangle->SetInputConnection(p_spline->GetOutputPort());
155  p_triangle->Update();
156 
157  vtkSmartPointer<vtkPolyData> p_cleaned = p_triangle->GetOutput();
158 
159  // Want flat ends on input and output network nodes. It is important to do this after smoothing
160  std::vector<boost::shared_ptr<VesselNode<DIM> > > nodes = mpNetwork->GetNodes();
161  c_vector<double, 3> box_axis = unit_vector<double>(3,0);
162 
163  for(unsigned idx=0; idx< nodes.size(); idx++)
164  {
165  if(nodes[idx]->GetFlowProperties()->IsInputNode() or nodes[idx]->GetFlowProperties()->IsOutputNode())
166  {
167  double radius = nodes[idx]->GetRadius()/mReferenceLength;
168  vtkSmartPointer<vtkBox> p_box = vtkSmartPointer<vtkBox>::New();
169  p_box->SetBounds(-1.1*radius, 0.0, -1.1*radius, 1.1*radius, - 1.1*radius, 1.1*radius);
170 
171  c_vector<double, 3> loc;
172  loc[0]= nodes[idx]->rGetLocation().GetLocation(mReferenceLength)[0];
173  loc[1]= nodes[idx]->rGetLocation().GetLocation(mReferenceLength)[1];
174  loc[2]=0.0;
175 
176  c_vector<double, 3> tangent;
177  tangent[0]= nodes[idx]->GetSegment(0)->GetOppositeNode(nodes[idx])->rGetLocation().GetLocation(mReferenceLength)[0] - loc[0];
178  tangent[1]= nodes[idx]->GetSegment(0)->GetOppositeNode(nodes[idx])->rGetLocation().GetLocation(mReferenceLength)[1] - loc[1];
179  tangent[2] = 0.0;
180  tangent /=norm_2(tangent);
181 
182  double rotation_angle = std::acos(inner_prod(box_axis, tangent))*(180.0/M_PI);
183  c_vector<double, 3> rotation_axis = VectorProduct(box_axis, tangent);
184 
185  vtkSmartPointer<vtkTransform> p_tranform = vtkSmartPointer<vtkTransform>::New();
186  if (std::abs(inner_prod(box_axis, tangent)) < 1.0 - 1.e-6)
187  {
188  p_tranform->RotateWXYZ(-rotation_angle, rotation_axis[0], rotation_axis[1], rotation_axis[2]);
189  }
190  else
191  {
192  p_tranform->RotateWXYZ(-rotation_angle, 0.0, 0.0, 1.0);
193  }
194  p_tranform->Translate(-loc[0],-loc[1], -loc[2]);
195  p_box->SetTransform(p_tranform);
196 
197  vtkSmartPointer<vtkClipPolyData> p_clipper = vtkSmartPointer<vtkClipPolyData>::New();
198  #if VTK_MAJOR_VERSION <= 5
199  p_clipper->SetInput(p_cleaned);
200  #else
201  p_clipper->SetInputData(p_cleaned);
202  #endif
203  p_clipper->SetClipFunction(p_box);
204  p_clipper->Update();
205 
206  // Assuming we have two points with connectivity 1, join them
207  std::vector<unsigned> edge_ids;
208  vtkSmartPointer<vtkIdList> p_cell_list = vtkSmartPointer<vtkIdList>::New();
209  for (unsigned idx =0; idx< p_clipper->GetOutput()->GetNumberOfPoints(); idx++)
210  {
211  p_clipper->GetOutput()->GetPointCells(idx, p_cell_list);
212  if(p_cell_list->GetNumberOfIds()==1)
213  {
214  edge_ids.push_back(idx);
215  }
216  }
217 
218  vtkSmartPointer<vtkLine> p_line = vtkSmartPointer<vtkLine>::New();
219  p_line->GetPointIds()->SetId(0, edge_ids[0]);
220  p_line->GetPointIds()->SetId(1, edge_ids[1]);
221  p_clipper->GetOutput()->GetLines()->InsertNextCell(p_line);
222 
223  p_cleaned->DeepCopy(p_clipper->GetOutput());
224  }
225  }
226  mpSurface = p_cleaned;
227  }
228  else
229  {
230  // Interpolate the network onto a regular grid
231  boost::shared_ptr<NetworkToImage<DIM> > p_net_to_image = NetworkToImage<DIM>::Create();
232  p_net_to_image->SetNetwork(mpNetwork);
233  p_net_to_image->SetGridSpacing(mSamplingGridSpacing);
234  p_net_to_image->SetPaddingFactors(0.1, 0.1, 0.1);
235  p_net_to_image->SetImageDimension(3);
236  p_net_to_image->Update();
237  mpImage = p_net_to_image->GetOutput();
238 
239  // Get the outer boundaries of the network
240  vtkSmartPointer<vtkMarchingCubes> p_cubes = vtkSmartPointer<vtkMarchingCubes>::New();
241  #if VTK_MAJOR_VERSION <= 5
242  p_cubes->SetInput(mpImage);
243  #else
244  p_cubes->SetInputData(mpImage);
245  #endif
246  p_cubes->SetValue(0, 1);
247 
248  vtkSmartPointer<vtkCleanPolyData> p_cleaner = vtkSmartPointer<vtkCleanPolyData>::New();
249  p_cleaner->SetInputConnection(p_cubes->GetOutputPort());
250 
251  vtkSmartPointer<vtkWindowedSincPolyDataFilter> p_smoother = vtkSmartPointer<vtkWindowedSincPolyDataFilter>::New();
252  p_smoother->SetInputConnection(p_cleaner->GetOutputPort());
253  p_smoother->SetNumberOfIterations(500.0);
254  p_smoother->BoundarySmoothingOn();
255  p_smoother->FeatureEdgeSmoothingOn();
256  p_smoother->SetFeatureAngle(30.0);
257  p_smoother->SetPassBand(0.03);
258  p_smoother->NonManifoldSmoothingOn();
259  p_smoother->NormalizeCoordinatesOn();
260  p_smoother->Update();
261 
262  vtkSmartPointer<vtkPolyData> p_cleaned = p_smoother->GetOutput();
263 
264  // Open any ends marked as inlet or outlet nodes
265  std::vector<boost::shared_ptr<VesselNode<DIM> > > nodes = mpNetwork->GetNodes();
266  c_vector<double, 3> box_axis = unit_vector<double>(3,0);
267  for(unsigned idx=0; idx< nodes.size(); idx++)
268  {
269  if(nodes[idx]->GetFlowProperties()->IsInputNode() or nodes[idx]->GetFlowProperties()->IsOutputNode())
270  {
271  double radius = nodes[idx]->GetRadius()/nodes[idx]->GetReferenceLengthScale();
272  vtkSmartPointer<vtkBox> p_box = vtkSmartPointer<vtkBox>::New();
273 
274  p_box->SetBounds(-1.1*radius, 0.0, -1.1*radius, 1.1*radius, - 1.1*radius, 1.1*radius);
275 
276  c_vector<double, 3> loc = nodes[idx]->rGetLocation().GetLocation(mReferenceLength);
277  c_vector<double, 3> tangent;
278  tangent = nodes[idx]->GetSegment(0)->GetOppositeNode(nodes[idx])->rGetLocation().GetLocation(mReferenceLength) - loc;
279  tangent /=norm_2(tangent);
280 
281  double rotation_angle = std::acos(inner_prod(box_axis, tangent))*(180.0/M_PI);
282  c_vector<double, 3> rotation_axis = VectorProduct(box_axis, tangent);
283 
284  vtkSmartPointer<vtkTransform> p_tranform = vtkSmartPointer<vtkTransform>::New();
285  if (std::abs(inner_prod(box_axis, tangent)) < 1.0 - 1.e-6)
286  {
287  p_tranform->RotateWXYZ(-rotation_angle, rotation_axis[0], rotation_axis[1], rotation_axis[2]);
288  }
289  else
290  {
291  p_tranform->RotateWXYZ(-rotation_angle, 0.0, 0.0, 1.0);
292  }
293  p_tranform->Translate(-loc[0],-loc[1], -loc[2]);
294  p_box->SetTransform(p_tranform);
295 
296  vtkSmartPointer<vtkClipPolyData> p_clipper = vtkSmartPointer<vtkClipPolyData>::New();
297  #if VTK_MAJOR_VERSION <= 5
298  p_clipper->SetInput(p_cleaned);
299  #else
300  p_clipper->SetInputData(p_cleaned);
301  #endif
302  p_clipper->SetClipFunction(p_box);
303  p_clipper->Update();
304  p_cleaned->DeepCopy(p_clipper->GetOutput());
305  }
306  }
307 
308  // Use vmtk to cap the surface
309 // vtkSmartPointer<vtkvmtkSimpleCapPolyData> p_capper = vtkSmartPointer<vtkvmtkSimpleCapPolyData>::New();
310 // p_capper->SetInputData(p_cleaned);
311 // p_capper->SetCellEntityIdsArrayName("CellEntityIds");
312 // p_capper->SetCellEntityIdOffset(1);
313 // p_capper->Update();
314 
315  vtkSmartPointer<vtkPolyDataNormals> p_normals = vtkSmartPointer<vtkPolyDataNormals>::New();
316 // p_normals->SetInputConnection(p_capper->GetOutputPort());
317  #if VTK_MAJOR_VERSION <= 5
318  p_normals->SetInput(p_cleaned);
319  #else
320  p_normals->SetInputData(p_cleaned);
321  #endif
322 
323  p_normals->AutoOrientNormalsOn();
324  p_normals->SplittingOff();
325  p_normals->ConsistencyOn();
326  p_normals->Update();
327 
328  mpSurface = p_normals->GetOutput();
329  }
330 }
331 
332 // Explicit instantiation
333 template class NetworkToSurface<2>;
334 template class NetworkToSurface<3>;
vtkSmartPointer< vtkImageData > mpImage
units::quantity< unit::length > mReferenceLength
units::quantity< unit::length > mSamplingGridSpacing
vtkSmartPointer< vtkImageData > GetSamplingImage()
units::quantity< unit::length > mSplineResamplingLength
vtkSmartPointer< vtkPolyData > mpSurface
void SetResamplingSplineSize(units::quantity< unit::length > splineResampleSize)
void SetResamplingGridSize(units::quantity< unit::length > sampleGridSize)
boost::shared_ptr< VesselNetwork< DIM > > mpNetwork
vtkSmartPointer< vtkPolyData > GetSurface()
static boost::shared_ptr< NetworkToSurface< DIM > > Create()
void SetVesselNetwork(boost::shared_ptr< VesselNetwork< DIM > > pNetwork)
static boost::shared_ptr< NetworkToImage< DIM > > Create()