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