37 #define _BACKWARD_BACKWARD_WARNING_H 1 //Cut out the vtk deprecated warning for now (gcc4.3) 40 #include <vtkPoints.h> 41 #include <vtkSmartPointer.h> 42 #include "RandomNumberGenerator.hpp" 43 #include "Polygon.hpp" 44 #include "BaseUnits.hpp" 46 #include "VoronoiGenerator.hpp" 58 template<
unsigned DIM>
60 mReferenceLength(
BaseUnits::Instance()->GetReferenceLengthScale())
64 template<
unsigned DIM>
70 template<
unsigned DIM>
77 EXCEPTION(
"This generator works in 3D only");
81 std::vector<units::quantity<unit::length> > extents = pPart->GetBoundingBox();
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;
100 for(
unsigned idx = 0; idx< numSeeds; idx++)
102 c_vector<double, DIM> location;
103 for(
unsigned jdx=0; jdx<DIM; jdx++)
105 location[jdx] = (dimensionless_extents[2*jdx] + (dimensionless_extents[2*jdx + 1] - dimensionless_extents[2*jdx])*RandomNumberGenerator::Instance()->ranf());
112 class tetgen::tetgenio mesher_input, mesher_output;
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++)
120 c_vector<double, DIM> seed_location = seeds[idx]->GetLocation(
mReferenceLength);
121 for (
unsigned jdx=0; jdx < DIM; jdx++)
123 mesher_input.pointlist[DIM*idx + jdx] = double(seed_location[jdx]);
127 tetgen::tetrahedralize((
char*)
"veeQ", &mesher_input, &mesher_output);
130 std::vector<boost::shared_ptr<Polygon<DIM> > > polygons;
131 for (
int n=0; n<mesher_output.numberofvedges; ++n)
133 tetgen::tetgenio::voroedge e = mesher_output.vedgelist[n];
136 double* u = &mesher_output.vpointlist[DIM*n0];
137 double* v = n1 == -1 ? e.vnormal : &mesher_output.vpointlist[DIM*n1];
142 double zero_location_tol = 20.0;
147 c_vector<double,DIM> intercept_1;
148 c_vector<double,DIM> intercept_2;
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);
165 p_polygon->AddVertex(p_vertex2);
166 polygons.push_back(p_polygon);
173 bool vert1_inside =
true;
174 bool vert2_inside =
true;
177 vert1_inside =
false;
181 vert1_inside =
false;
185 vert2_inside =
false;
189 vert2_inside =
false;
195 vert1_inside =
false;
199 vert2_inside =
false;
203 if(vert1_inside && !vert2_inside)
206 c_vector<double,DIM> gap = intercept_1-p_vertex2->GetLocation(
mReferenceLength);
209 else if(vert2_inside && !vert1_inside)
212 c_vector<double,DIM> gap = intercept_1-p_vertex1->GetLocation(
mReferenceLength);
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)
224 if(norm_2(tv2) > zero_location_tol)
236 p_polygon->AddVertex(p_vertex2);
237 polygons.push_back(p_polygon);
245 for (
unsigned idx=0; idx<polygons.size(); idx++)
247 p_tesselation->AddPolygon(polygons[idx]);
250 return p_tesselation;
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)
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()
units::quantity< unit::length > mReferenceLength