Chaste  Build::
BoundaryExtractor.cpp
1 /*
2  Copyright (c) 2005-2015, University of Oxford.
3  All rights reserved.
4 
5  University of Oxford means the Chancellor, Masters and Scholars of the
6  University of Oxford, having an administrative office at Wellington
7  Square, Oxford OX1 2JD, UK.
8 
9  This file is part of Chaste.
10 
11  Redistribution and use in source and binary forms, with or without
12  modification, are permitted provided that the following conditions are met:
13  * Redistributions of source code must retain the above copyright notice,
14  this list of conditions and the following disclaimer.
15  * Redistributions in binary form must reproduce the above copyright notice,
16  this list of conditions and the following disclaimer in the documentation
17  and/or other materials provided with the distribution.
18  * Neither the name of the University of Oxford nor the names of its
19  contributors may be used to endorse or promote products derived from this
20  software without specific prior written permission.
21 
22  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
23  AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
24  IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
25  ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
26  LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
27  CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
28  GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
29  HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
30  LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
31  OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
32 
33  */
34 
35 #include <vtkFeatureEdges.h>
36 #include <vtkCleanPolyData.h>
37 #include <vtkTriangleFilter.h>
38 #include <vtkStripper.h>
39 #include <vtkSplineFilter.h>
40 #include <vtkConnectivityFilter.h>
41 #include <vtkGeometryFilter.h>
42 #include "Exception.hpp"
43 #include "BoundaryExtractor.hpp"
44 
46  : mpInputSurface(),
47  mpOutputSurface(),
48  mSmoothingLength(1),
49  mDoSmoothing(true)
50 {
51 
52 }
53 
54 boost::shared_ptr<BoundaryExtractor> BoundaryExtractor::Create()
55 {
56  MAKE_PTR(BoundaryExtractor, pSelf);
57  return pSelf;
58 }
59 
61 {
62 
63 }
64 
65 vtkSmartPointer<vtkPolyData> BoundaryExtractor::GetOutput()
66 {
67  if(mpOutputSurface)
68  {
69  return(mpOutputSurface);
70  }
71  else
72  {
73  EXCEPTION("No output set. Did you run 'Update()' ?");
74  }
75 }
76 
77 void BoundaryExtractor::SetInput(vtkSmartPointer<vtkPolyData> pInput)
78 {
79  mpInputSurface = pInput;
80 }
81 
83 {
84  mSmoothingLength = value;
85 }
86 
87 void BoundaryExtractor::SetDoSmoothing(bool doSmoothing)
88 {
89  mDoSmoothing = doSmoothing;
90 }
91 
93 {
94  if(!mpInputSurface)
95  {
96  EXCEPTION("No input set.");
97  }
98 
99  vtkSmartPointer<vtkFeatureEdges> p_features = vtkSmartPointer<vtkFeatureEdges>::New();
100  #if VTK_MAJOR_VERSION <= 5
101  p_features->SetInput(mpInputSurface);
102  #else
103  p_features->SetInputData(mpInputSurface);
104  #endif
105 
106  vtkSmartPointer<vtkCleanPolyData> p_clean = vtkSmartPointer<vtkCleanPolyData>::New();
107  p_clean->SetInputConnection(p_features->GetOutputPort());
108 
109  vtkSmartPointer<vtkTriangleFilter> p_triangle = vtkSmartPointer<vtkTriangleFilter>::New();
110  p_triangle->SetInputConnection(p_clean->GetOutputPort());
111  p_triangle->Update();
112 
113  if(mDoSmoothing)
114  {
115  vtkSmartPointer<vtkStripper> p_stripper = vtkSmartPointer<vtkStripper>::New();
116  p_stripper->SetInputConnection(p_triangle->GetOutputPort());
117 
118  vtkSmartPointer<vtkSplineFilter> p_spline = vtkSmartPointer<vtkSplineFilter>::New();
119  p_spline->SetInputConnection(p_stripper->GetOutputPort());
120  p_spline->SetSubdivideToLength();
121  p_spline->SetLength(mSmoothingLength);
122  p_spline->Update();
123 
124  vtkSmartPointer<vtkGeometryFilter> p_geom = vtkSmartPointer<vtkGeometryFilter>::New();
125  p_geom->SetInputConnection(p_spline->GetOutputPort());
126  p_geom->Update();
127 
128  vtkSmartPointer<vtkTriangleFilter> p_triangle2 = vtkSmartPointer<vtkTriangleFilter>::New();
129  p_triangle2->SetInputConnection(p_geom->GetOutputPort());
130  p_triangle2->Update();
131 
132  vtkSmartPointer<vtkCleanPolyData> p_clean2 = vtkSmartPointer<vtkCleanPolyData>::New();
133  p_clean2->SetInputConnection(p_triangle2->GetOutputPort());
134  p_clean2->Update();
135 
136  mpOutputSurface = p_clean2->GetOutput();
137  }
138  else
139  {
140  mpOutputSurface = p_triangle->GetOutput();
141  }
142 }
void SetSmoothingLength(double value)
vtkSmartPointer< vtkPolyData > mpInputSurface
vtkSmartPointer< vtkPolyData > GetOutput()
vtkSmartPointer< vtkPolyData > mpOutputSurface
static boost::shared_ptr< BoundaryExtractor > Create()
void SetDoSmoothing(bool doSmoothing)
void SetInput(vtkSmartPointer< vtkPolyData > pInputSurface)