Chaste  Build::
Polygon.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 vtk deprecated warning
37 #include <vtkIdTypeArray.h>
38 #include <vtkTriangle.h>
39 #include <vtkDoubleArray.h>
40 #include "Exception.hpp"
41 #include "Polygon.hpp"
42 
43 template<unsigned DIM>
44 Polygon<DIM>::Polygon(std::vector<boost::shared_ptr<DimensionalChastePoint<DIM> > > vertices) :
45  mVertices(vertices),
46  mReferenceLength(BaseUnits::Instance()->GetReferenceLengthScale())
47 {
48 
49 }
50 
51 template<unsigned DIM>
52 Polygon<DIM>::Polygon(boost::shared_ptr<DimensionalChastePoint<DIM> > pVertex) :
53  mVertices(),
54  mReferenceLength(BaseUnits::Instance()->GetReferenceLengthScale())
55 {
56  mVertices.push_back(pVertex);
57 }
58 
59 template<unsigned DIM>
60 boost::shared_ptr<Polygon<DIM> > Polygon<DIM>::Create(std::vector<boost::shared_ptr<DimensionalChastePoint<DIM> > > vertices)
61 {
62  MAKE_PTR_ARGS(Polygon<DIM> , pSelf, (vertices));
63  return pSelf;
64 }
65 
66 template<unsigned DIM>
67 boost::shared_ptr<Polygon<DIM> > Polygon<DIM>::Create(boost::shared_ptr<DimensionalChastePoint<DIM> > pVertex)
68 {
69  MAKE_PTR_ARGS(Polygon<DIM> , pSelf, (pVertex));
70  return pSelf;
71 }
72 
73 template<unsigned DIM>
75 {
76 }
77 
78 template<unsigned DIM>
79 void Polygon<DIM>::AddVertices(std::vector<boost::shared_ptr<DimensionalChastePoint<DIM> > > vertices)
80 {
81  // Add the new vertices
82  mVertices.insert(mVertices.end(), vertices.begin(), vertices.end());
83 }
84 
85 template<unsigned DIM>
86 void Polygon<DIM>::AddVertex(boost::shared_ptr<DimensionalChastePoint<DIM> > pVertex)
87 {
88  mVertices.push_back(pVertex);
89 }
90 
91 template<unsigned DIM>
92 boost::shared_ptr<DimensionalChastePoint<DIM> > Polygon<DIM>::GetVertex(unsigned idx)
93 {
94  if(idx >= mVertices.size())
95  {
96  EXCEPTION("Requested vertex index out of range");
97  }
98  else
99  {
100  return mVertices[idx];
101  }
102 }
103 
104 template<unsigned DIM>
105 std::vector<boost::shared_ptr<DimensionalChastePoint<DIM> > > Polygon<DIM>::GetVertices()
106 {
107  return mVertices;
108 }
109 
110 template<unsigned DIM>
111 std::pair<vtkSmartPointer<vtkPoints>, vtkSmartPointer<vtkIdTypeArray> > Polygon<DIM>::GetVtkVertices()
112 {
113  vtkSmartPointer<vtkIdTypeArray> p_vertexIds = vtkSmartPointer<vtkIdTypeArray>::New();
114  vtkSmartPointer<vtkPoints> p_vertices = vtkSmartPointer<vtkPoints>::New();
115  p_vertices->SetNumberOfPoints(mVertices.size());
116  p_vertexIds->SetNumberOfTuples(mVertices.size());
117  for (vtkIdType idx = 0; idx < vtkIdType(mVertices.size()); idx++)
118  {
119  c_vector<double, DIM> vertex_location = mVertices[idx]->GetLocation(mReferenceLength);
120  if(DIM==3)
121  {
122  p_vertices->SetPoint(idx, vertex_location[0], vertex_location[1], vertex_location[2]);
123  }
124  else
125  {
126  p_vertices->SetPoint(idx, vertex_location[0], vertex_location[1], 0.0);
127  }
128  p_vertexIds->SetValue(idx, idx);
129  }
130  return std::pair<vtkSmartPointer<vtkPoints>, vtkSmartPointer<vtkIdTypeArray> >(p_vertices, p_vertexIds);
131 }
132 
133 template<unsigned DIM>
134 vtkSmartPointer<vtkPolygon> Polygon<DIM>::GetVtkPolygon()
135 {
136  vtkSmartPointer<vtkPoints> p_vertices = GetVtkVertices().first;
137  vtkSmartPointer<vtkPolygon> p_polygon = vtkSmartPointer<vtkPolygon>::New();
138  p_polygon->GetPoints()->SetNumberOfPoints(p_vertices->GetNumberOfPoints());
139  for (vtkIdType idx = 0; idx < p_vertices->GetNumberOfPoints(); idx++)
140  {
141  p_polygon->GetPoints()->SetPoint(idx, p_vertices->GetPoint(idx));
142  }
143  return p_polygon;
144 }
145 
146 template<unsigned DIM>
148 {
149  c_vector<double, 3> centroid;
150  std::pair<vtkSmartPointer<vtkPoints>, vtkSmartPointer<vtkIdTypeArray> > vertex_data = GetVtkVertices();
151  vtkPolygon::ComputeCentroid(vertex_data.second, vertex_data.first, &centroid[0]);
152 
153  if(DIM==3)
154  {
156  }
157  else
158  {
159  return DimensionalChastePoint<DIM>(centroid[0], centroid[1], 0.0, mReferenceLength);
160  }
161 
162 }
163 
164 template<unsigned DIM>
165 std::vector<units::quantity<unit::length> > Polygon<DIM>::GetBoundingBox()
166 {
167  vtkSmartPointer<vtkPolygon> p_polygon = GetVtkPolygon();
168  c_vector<double, 6> box;
169  p_polygon->GetPoints()->GetBounds(&box[0]);
170 
171  std::vector<units::quantity<unit::length> > box_vector(6);
172  for(unsigned idx=0; idx<6; idx++)
173  {
174  box_vector[idx] = box[idx] * mReferenceLength;
175  }
176  return box_vector;
177 }
178 
179 template<unsigned DIM>
180 units::quantity<unit::length> Polygon<DIM>::GetDistance(const DimensionalChastePoint<DIM>& location)
181 {
182  double point[3];
183  for (unsigned idx = 0; idx < DIM; idx++)
184  {
185  point[idx] = location.GetLocation(mReferenceLength)[idx];
186  }
187  if(DIM==2)
188  {
189  point[2] = 0.0;
190  }
191 
192  vtkSmartPointer<vtkPlane> p_plane = GetPlane();
193  double distance = p_plane->DistanceToPlane(&point[0]);
194  return distance * mReferenceLength;
195 }
196 
197 template<unsigned DIM>
198 units::quantity<unit::length> Polygon<DIM>::GetDistanceToEdges(const DimensionalChastePoint<DIM>& location)
199 {
200  double point[3];
201  for (unsigned idx = 0; idx < DIM; idx++)
202  {
203  point[idx] = location.GetLocation(mReferenceLength)[idx];
204  }
205  if(DIM==2)
206  {
207  point[2] = 0.0;
208  }
209 
210  vtkSmartPointer<vtkPolygon> p_polygon = GetVtkPolygon();
211 
212  double bounds[6];
213  double closest[6];
214  p_polygon->GetPoints()->GetBounds(bounds);
215 
216  double distance = p_polygon->DistanceToPolygon(
217  point, p_polygon->GetPoints()->GetNumberOfPoints(),
218  static_cast<double*>(p_polygon->GetPoints()->GetData()->GetVoidPointer(0)), bounds, closest);
219 
220  return distance * mReferenceLength;
221 }
222 
223 template<unsigned DIM>
224 vtkSmartPointer<vtkPlane> Polygon<DIM>::GetPlane()
225 {
226  vtkSmartPointer<vtkPlane> p_plane = vtkSmartPointer<vtkPlane>::New();
227  c_vector<double, DIM> centroid = GetCentroid().GetLocation(mReferenceLength);
228  c_vector<double, DIM> normal = GetNormal();
229  if(DIM==3)
230  {
231  p_plane->SetOrigin(centroid[0], centroid[1], centroid[2]);
232  p_plane->SetNormal(normal[0], normal[1], normal[2]);
233  }
234  else
235  {
236  p_plane->SetOrigin(centroid[0], centroid[1], 0.0);
237  p_plane->SetNormal(normal[0], normal[1], 0.0);
238  }
239  return p_plane;
240 }
241 
242 template<unsigned DIM>
243 c_vector<double, DIM> Polygon<DIM>::GetNormal()
244 {
245  if (mVertices.size() < 3)
246  {
247  EXCEPTION("At least 3 vertices are required to generate a normal.");
248  }
249 
250  std::pair<vtkSmartPointer<vtkPoints>, vtkSmartPointer<vtkIdTypeArray> > vertex_data = GetVtkVertices();
251  double loc1[3];
252  double loc2[3];
253  double loc3[3];
254  vertex_data.first->GetPoint(0, loc1);
255  vertex_data.first->GetPoint(1, loc2);
256  vertex_data.first->GetPoint(2, loc3);
257  c_vector<double, 3> in_normal;
258  vtkTriangle::ComputeNormal(loc1, loc2, loc3, &in_normal[0]);
259  if(DIM==3)
260  {
261  return in_normal;
262  }
263  else
264  {
265  c_vector<double, 2> normal;
266  normal[0] = in_normal[0];
267  normal[1] = in_normal[1];
268  return normal;
269  }
270 }
271 
272 template<unsigned DIM>
274 {
275  bool contains_point = false;
276  if (mVertices.size() >= 3)
277  {
278  c_vector<double, DIM> vertex_location = location.GetLocation(mReferenceLength);
279  double point[3];
280  for (unsigned idx = 0; idx < DIM; idx++)
281  {
282  point[idx] = vertex_location[idx];
283  }
284  if(DIM==2)
285  {
286  point[2] = 0.0;
287  }
288 
289  vtkSmartPointer<vtkPolygon> p_polygon = GetVtkPolygon();
290  double n[3];
291  p_polygon->ComputeNormal(p_polygon->GetPoints()->GetNumberOfPoints(),
292  static_cast<double*>(p_polygon->GetPoints()->GetData()->GetVoidPointer(0)), n);
293 
294  double bounds[6];
295  p_polygon->GetPoints()->GetBounds(bounds);
296 
297  int contains = p_polygon->PointInPolygon(
298  point, p_polygon->GetPoints()->GetNumberOfPoints(),
299  static_cast<double*>(p_polygon->GetPoints()->GetData()->GetVoidPointer(0)), bounds, n);
300 
301  if (contains == 1)
302  {
303  contains_point = true;
304  }
305  }
306  return contains_point;
307 }
308 
309 template<unsigned DIM>
310 void Polygon<DIM>::ReplaceVertex(unsigned idx, boost::shared_ptr<DimensionalChastePoint<DIM> > pVertex)
311 {
312  if(idx >= mVertices.size())
313  {
314  EXCEPTION("Requested vertex index out of range");
315  }
316  else
317  {
318  mVertices[idx] = pVertex;
319  }
320 }
321 
322 template<unsigned DIM>
323 void Polygon<DIM>::RotateAboutAxis(c_vector<double, 3> axis, double angle)
324 {
325  for(unsigned idx=0; idx<mVertices.size(); idx++)
326  {
327  mVertices[idx]->RotateAboutAxis(axis, angle);
328  }
329 }
330 
331 template<unsigned DIM>
333 {
334  for(unsigned idx=0; idx<mVertices.size(); idx++)
335  {
336  mVertices[idx]->Translate(translationVector);
337  }
338 }
339 
340 // Explicit instantiation
341 template class Polygon<2>;
342 template class Polygon<3>;
343 
Polygon(std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > > vertices)
Definition: Polygon.cpp:44
vtkSmartPointer< vtkPolygon > GetVtkPolygon()
Definition: Polygon.cpp:134
~Polygon()
Definition: Polygon.cpp:74
bool ContainsPoint(const DimensionalChastePoint< DIM > &rLocation)
Definition: Polygon.cpp:273
void AddVertex(boost::shared_ptr< DimensionalChastePoint< DIM > > pVertex)
Definition: Polygon.cpp:86
void Translate(DimensionalChastePoint< DIM > translationVector)
Definition: Polygon.cpp:332
units::quantity< unit::length > GetDistance(const DimensionalChastePoint< DIM > &rLocation)
Definition: Polygon.cpp:180
static boost::shared_ptr< Polygon > Create(std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > > vertices)
Definition: Polygon.cpp:60
void ReplaceVertex(unsigned idx, boost::shared_ptr< DimensionalChastePoint< DIM > > pVertex)
Definition: Polygon.cpp:310
std::pair< vtkSmartPointer< vtkPoints >, vtkSmartPointer< vtkIdTypeArray > > GetVtkVertices()
Definition: Polygon.cpp:111
c_vector< double, DIM > GetLocation(units::quantity< unit::length > scale)
std::vector< units::quantity< unit::length > > GetBoundingBox()
Definition: Polygon.cpp:165
std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > > mVertices
Definition: Polygon.hpp:60
void AddVertices(std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > > vertices)
Definition: Polygon.cpp:79
units::quantity< unit::length > mReferenceLength
Definition: Polygon.hpp:65
DimensionalChastePoint< DIM > GetCentroid()
Definition: Polygon.cpp:147
units::quantity< unit::length > GetDistanceToEdges(const DimensionalChastePoint< DIM > &rLocation)
Definition: Polygon.cpp:198
vtkSmartPointer< vtkPlane > GetPlane()
Definition: Polygon.cpp:224
boost::shared_ptr< DimensionalChastePoint< DIM > > GetVertex(unsigned idx)
Definition: Polygon.cpp:92
void RotateAboutAxis(c_vector< double, 3 > axis, double angle)
Definition: Polygon.cpp:323
std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > > GetVertices()
Definition: Polygon.cpp:105
c_vector< double, DIM > GetNormal()
Definition: Polygon.cpp:243