Chaste  Build::
VoronoiGenerator.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 <math.h>
37 #define _BACKWARD_BACKWARD_WARNING_H 1 //Cut out the vtk deprecated warning for now (gcc4.3)
38 #include <vtkBox.h>
39 #include <vtkTetra.h>
40 #include <vtkPoints.h>
41 #include <vtkSmartPointer.h>
42 #include "RandomNumberGenerator.hpp"
43 #include "Polygon.hpp"
44 #include "BaseUnits.hpp"
45 
46 #include "VoronoiGenerator.hpp"
47 
48 // Hang Si's tetgen.
49 
50 #define REAL double
51 #define VOID void
52 #include "tetgen.h"
53 #undef REAL
54 #undef VOID
55 
56 struct triangulateio;
57 
58 template<unsigned DIM>
60  mReferenceLength(BaseUnits::Instance()->GetReferenceLengthScale())
61 {
62 }
63 
64 template<unsigned DIM>
66 {
67 
68 }
69 
70 template<unsigned DIM>
71 boost::shared_ptr<Part<DIM> > VoronoiGenerator<DIM>::Generate(boost::shared_ptr<Part<DIM> > pPart,
72  std::vector<boost::shared_ptr<DimensionalChastePoint<DIM> > > seeds,
73  unsigned numSeeds)
74 {
75  if(DIM==2)
76  {
77  EXCEPTION("This generator works in 3D only");
78  }
79 
80  boost::shared_ptr<Part<DIM> > p_tesselation = Part<DIM>::Create();
81  std::vector<units::quantity<unit::length> > extents = pPart->GetBoundingBox();
82  double x_min = double(extents[0]/ mReferenceLength);
83  double x_max = double(extents[1]/ mReferenceLength);
84  double y_min = double(extents[2]/ mReferenceLength);
85  double y_max = double(extents[3]/ mReferenceLength);
86  double z_min = double(extents[4]/ mReferenceLength);
87  double z_max = double(extents[5]/ mReferenceLength);
88 
89  c_vector<double, 6> dimensionless_extents;
90  dimensionless_extents[0] = x_min;
91  dimensionless_extents[1] = x_max;
92  dimensionless_extents[2] = y_min;
93  dimensionless_extents[3] = y_max;
94  dimensionless_extents[4] = z_min;
95  dimensionless_extents[5] = z_max;
96 
97  // If no seeds have been provided generate some random ones
98  if(seeds.size() == 0)
99  {
100  for(unsigned idx = 0; idx< numSeeds; idx++)
101  {
102  c_vector<double, DIM> location;
103  for(unsigned jdx=0; jdx<DIM; jdx++)
104  {
105  location[jdx] = (dimensionless_extents[2*jdx] + (dimensionless_extents[2*jdx + 1] - dimensionless_extents[2*jdx])*RandomNumberGenerator::Instance()->ranf());
106  }
107  seeds.push_back(DimensionalChastePoint<DIM>::Create(location, mReferenceLength));
108  }
109  }
110 
111  // Use tetgen to make the tessellation
112  class tetgen::tetgenio mesher_input, mesher_output;
113 
114  mesher_input.initialize();
115  mesher_output.initialize();
116  mesher_input.numberofpoints = seeds.size();
117  mesher_input.pointlist = new double[mesher_input.numberofpoints * DIM];
118  for (unsigned idx=0; idx<seeds.size(); idx++)
119  {
120  c_vector<double, DIM> seed_location = seeds[idx]->GetLocation(mReferenceLength);
121  for (unsigned jdx=0; jdx < DIM; jdx++)
122  {
123  mesher_input.pointlist[DIM*idx + jdx] = double(seed_location[jdx]);
124  }
125  }
126 
127  tetgen::tetrahedralize((char*)"veeQ", &mesher_input, &mesher_output);
128 
129  // Create 2-point polygons corresponding to each edge
130  std::vector<boost::shared_ptr<Polygon<DIM> > > polygons;
131  for (int n=0; n<mesher_output.numberofvedges; ++n)
132  {
133  tetgen::tetgenio::voroedge e = mesher_output.vedgelist[n];
134  int n0 = e.v1;
135  int n1 = e.v2;
136  double* u = &mesher_output.vpointlist[DIM*n0];
137  double* v = n1 == -1 ? e.vnormal : &mesher_output.vpointlist[DIM*n1]; // -1 indicates ray
138  boost::shared_ptr<DimensionalChastePoint<DIM> > p_vertex1 = DimensionalChastePoint<DIM>::Create(u[0], u[1], u[2], mReferenceLength);
139  boost::shared_ptr<DimensionalChastePoint<DIM> > p_vertex2 = DimensionalChastePoint<DIM>::Create(v[0], v[1], v[2], mReferenceLength);
140 
141  // Only include edges that are inside the domain - use the bounding box.
142  double zero_location_tol = 20.0;
143  double t1;
144  double t2;
145  int plane1;
146  int plane2;
147  c_vector<double,DIM> intercept_1;
148  c_vector<double,DIM> intercept_2;
149 
150  int in_box = vtkBox::IntersectWithLine(&dimensionless_extents[0], &p_vertex1->GetLocation(mReferenceLength)[0],
151  &p_vertex2->GetLocation(mReferenceLength)[0], t1, t2, &intercept_1[0], &intercept_2[0], plane1, plane2);
152 
153  // If the line is not outside the box
154  if(in_box!=0)
155  {
156  // if the line does not intersect the box it must be fully inside
157  if(in_box==-1)
158  {
159  // Check that we haven't hit the zero location and remove very short edges
160  if((norm_2(p_vertex1->GetLocation(mReferenceLength)) > zero_location_tol) && (norm_2(p_vertex2->GetLocation(mReferenceLength)) > zero_location_tol))
161  {
162  if(norm_2(p_vertex1->GetLocation(mReferenceLength) - p_vertex2->GetLocation(mReferenceLength)) > 2.0 * zero_location_tol)
163  {
164  boost::shared_ptr<Polygon<DIM> > p_polygon = Polygon<DIM>::Create(p_vertex1);
165  p_polygon->AddVertex(p_vertex2);
166  polygons.push_back(p_polygon);
167  }
168  }
169  }
170  else
171  {
172  // Find where the points are
173  bool vert1_inside = true;
174  bool vert2_inside = true;
175  if(p_vertex1->GetLocation(mReferenceLength)[0] < dimensionless_extents[0] || p_vertex1->GetLocation(mReferenceLength)[0] > dimensionless_extents[1])
176  {
177  vert1_inside = false;
178  }
179  if(p_vertex1->GetLocation(mReferenceLength)[1] < dimensionless_extents[2] || p_vertex1->GetLocation(mReferenceLength)[1] > dimensionless_extents[3])
180  {
181  vert1_inside = false;
182  }
183  if(p_vertex2->GetLocation(mReferenceLength)[0] < dimensionless_extents[0] || p_vertex2->GetLocation(mReferenceLength)[0] > dimensionless_extents[1])
184  {
185  vert2_inside = false;
186  }
187  if(p_vertex2->GetLocation(mReferenceLength)[1] < dimensionless_extents[2] || p_vertex2->GetLocation(mReferenceLength)[1] > dimensionless_extents[3])
188  {
189  vert2_inside = false;
190  }
191  if(DIM==3)
192  {
193  if(p_vertex1->GetLocation(mReferenceLength)[2] < dimensionless_extents[4] || p_vertex1->GetLocation(mReferenceLength)[2] > dimensionless_extents[5])
194  {
195  vert1_inside = false;
196  }
197  if(p_vertex2->GetLocation(mReferenceLength)[2] < dimensionless_extents[4] || p_vertex2->GetLocation(mReferenceLength)[2] > dimensionless_extents[5])
198  {
199  vert2_inside = false;
200  }
201  }
202 
203  if(vert1_inside && !vert2_inside)
204  {
205  // move vert 2 to intsersection point
206  c_vector<double,DIM> gap = intercept_1-p_vertex2->GetLocation(mReferenceLength);
207  p_vertex2->Translate(DimensionalChastePoint<DIM>(gap, mReferenceLength));
208  }
209  else if(vert2_inside && !vert1_inside)
210  {
211  // move vert 1 to intersection point
212  c_vector<double,DIM> gap = intercept_1-p_vertex1->GetLocation(mReferenceLength);
213  p_vertex1->Translate(DimensionalChastePoint<DIM>(gap, mReferenceLength));
214  }
215  else
216  {
217  // Otherwise both points are either on or outside
218  c_vector<double,DIM> tv1 = intercept_1-p_vertex1->GetLocation(mReferenceLength);
219  c_vector<double,DIM> tv2 = intercept_2-p_vertex2->GetLocation(mReferenceLength);
220  if(norm_2(tv1) > zero_location_tol)
221  {
222  p_vertex1->Translate(DimensionalChastePoint<DIM>(tv1, mReferenceLength));
223  }
224  if(norm_2(tv2) > zero_location_tol)
225  {
226  p_vertex2->Translate(DimensionalChastePoint<DIM>(tv2, mReferenceLength));
227  }
228  }
229 
230  // Check that we haven't hit the zero location and remove very short edges
231  if((norm_2(p_vertex1->GetLocation(mReferenceLength)) > zero_location_tol) && (norm_2(p_vertex2->GetLocation(mReferenceLength)) > zero_location_tol))
232  {
233  if(norm_2(p_vertex1->GetLocation(mReferenceLength) - p_vertex2->GetLocation(mReferenceLength)) > 2.0 * zero_location_tol)
234  {
235  boost::shared_ptr<Polygon<DIM> > p_polygon = Polygon<DIM>::Create(p_vertex1);
236  p_polygon->AddVertex(p_vertex2);
237  polygons.push_back(p_polygon);
238  }
239  }
240  }
241  }
242  }
243 
244  // Add the polygons to the part
245  for (unsigned idx=0; idx<polygons.size(); idx++)
246  {
247  p_tesselation->AddPolygon(polygons[idx]);
248  }
249 
250  return p_tesselation;
251 }
252 
253 //Explicit instantiation
254 template class VoronoiGenerator<2> ;
255 template class VoronoiGenerator<3> ;
boost::shared_ptr< Part< DIM > > Generate(boost::shared_ptr< Part< DIM > > pPart, std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > > seeds=std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > >(), unsigned numSeeds=100)
static boost::shared_ptr< Polygon > Create(std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > > vertices)
Definition: Polygon.cpp:60
static boost::shared_ptr< DimensionalChastePoint< DIM > > Create(double x, double y, double z, units::quantity< unit::length > referenceLength)
static boost::shared_ptr< Part< DIM > > Create()
Definition: Part.cpp:63
units::quantity< unit::length > mReferenceLength
Definition: Part.hpp:60