36 #include <boost/lexical_cast.hpp> 37 #define _BACKWARD_BACKWARD_WARNING_H 1 //Cut out the vtk deprecated warning 38 #include <vtkIdList.h> 39 #include <vtkCellArray.h> 41 #include <vtkPoints.h> 42 #include "Exception.hpp" 43 #include "Warnings.hpp" 45 #include "Polygon.hpp" 46 #include "Element.hpp" 47 #include "UblasVectorInclude.hpp" 48 #include "AbstractTetrahedralMesh.hpp" 49 #include "VtkMeshWriter.hpp" 50 #include "DiscreteContinuumMesh.hpp" 51 #include "DiscreteContinuumMeshGenerator.hpp" 52 #include "BaseUnits.hpp" 54 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
69 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
76 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
82 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
87 EXCEPTION(
"No mesh has been generated.");
92 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
98 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
104 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
110 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
116 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
122 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
128 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
135 if (ELEMENT_DIM == 2)
145 std::vector<units::quantity<unit::length> > bounding_box =
mpDomain->GetBoundingBox();
146 if (units::abs(bounding_box[4]) < 1.e-12*unit::metres && units::abs(bounding_box[5]) < 1.e-12*unit::metres)
152 EXCEPTION(
"2D meshing is only supported for parts with z=0.");
161 EXCEPTION(
"2d meshing is not supported for STL files.");
170 std::vector<units::quantity<unit::length> > bounding_box =
mpDomain->GetBoundingBox();
171 if (units::abs(bounding_box[4]) < 1.e-12*unit::metres && units::abs(bounding_box[5]) < 1.e-12*unit::metres)
173 EXCEPTION(
"The part is two-dimensional, use the 2D meshing functionality.");
191 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
194 struct triangulateio mesher_input, mesher_output;
201 std::vector<DimensionalChastePoint<SPACE_DIM> > vertex_locations =
mpDomain->GetVertexLocations();
202 unsigned num_vertices = vertex_locations.size();
203 mesher_input.pointlist = (
double *) malloc(num_vertices * 2 *
sizeof(
double));
204 mesher_input.numberofpoints = int(num_vertices);
205 for (
unsigned idx = 0; idx < num_vertices; idx++)
207 c_vector<double, SPACE_DIM> vertex_location = vertex_locations[idx].GetLocation(
mReferenceLength);
208 for (
unsigned jdx = 0; jdx < 2; jdx++)
210 mesher_input.pointlist[2 * idx + jdx] = vertex_location[jdx];
214 std::vector<std::pair<unsigned, unsigned> > segments =
mpDomain->GetSegmentIndices();
215 unsigned num_segments = segments.size();
216 mesher_input.segmentlist = (
int *) malloc(num_segments * 2 *
sizeof(
int));
217 mesher_input.numberofsegments = int(num_segments);
218 for (
unsigned idx = 0; idx < num_segments; idx++)
220 mesher_input.segmentlist[2 * idx] = int(segments[idx].first);
221 mesher_input.segmentlist[2 * idx + 1] = int(segments[idx].second);
226 unsigned num_points = mpVtkDomain->GetNumberOfPoints();
227 mesher_input.pointlist = (
double *) malloc(num_points * 2 *
sizeof(
double));
228 mesher_input.numberofpoints = int(num_points);
229 for (
unsigned idx = 0; idx < num_points; idx++)
231 for (
unsigned jdx = 0; jdx < 2; jdx++)
233 mesher_input.pointlist[2 * idx + jdx] = mpVtkDomain->GetPoints()->GetPoint(idx)[jdx];
237 unsigned num_segments = mpVtkDomain->GetNumberOfLines();
238 mpVtkDomain->GetLines()->InitTraversal();
239 mesher_input.segmentlist = (
int *) malloc(num_segments * 2 *
sizeof(
int));
240 mesher_input.numberofsegments = int(num_segments);
241 vtkSmartPointer<vtkIdList> p_id_list = vtkSmartPointer<vtkIdList>::New();
242 for (
unsigned idx = 0; idx < num_segments; idx++)
244 mpVtkDomain->GetLines()->GetNextCell(p_id_list);
245 mesher_input.segmentlist[2 * idx] = int(p_id_list->GetId(0));
246 mesher_input.segmentlist[2 * idx + 1] = int(p_id_list->GetId(1));
251 unsigned num_points = mpVtkDomain->GetNumberOfPoints();
252 std::vector<DimensionalChastePoint<SPACE_DIM> > vertex_locations =
mpDomain->GetVertexLocations();
253 unsigned num_vertices = vertex_locations.size();
255 mesher_input.pointlist = (
double *) malloc((num_points+num_vertices) * 2 *
sizeof(double));
256 mesher_input.numberofpoints = int(num_points+num_vertices);
257 for (
unsigned idx = 0; idx < num_points; idx++)
259 for (
unsigned jdx = 0; jdx < 2; jdx++)
261 mesher_input.pointlist[2 * idx + jdx] = mpVtkDomain->GetPoints()->GetPoint(idx)[jdx];
264 for (
unsigned idx = 0; idx < num_vertices; idx++)
266 c_vector<double, SPACE_DIM> vertex_location = vertex_locations[idx].GetLocation(
mReferenceLength);
267 for (
unsigned jdx = 0; jdx < 2; jdx++)
269 mesher_input.pointlist[2 * idx + jdx + 2*num_points] = vertex_location[jdx];
273 unsigned num_segments = mpVtkDomain->GetNumberOfLines();
274 std::vector<std::pair<unsigned, unsigned> > segments =
mpDomain->GetSegmentIndices();
275 unsigned num_part_segments = segments.size();
276 mpVtkDomain->GetLines()->InitTraversal();
278 mesher_input.segmentlist = (
int *) malloc((num_segments+num_part_segments) * 2 *
sizeof(int));
279 mesher_input.numberofsegments = int(num_segments+num_part_segments);
280 vtkSmartPointer<vtkIdList> p_id_list = vtkSmartPointer<vtkIdList>::New();
281 for (
unsigned idx = 0; idx < num_segments; idx++)
283 mpVtkDomain->GetLines()->GetNextCell(p_id_list);
284 mesher_input.segmentlist[2 * idx] = int(p_id_list->GetId(0));
285 mesher_input.segmentlist[2 * idx + 1] = int(p_id_list->GetId(1));
287 for (
unsigned idx = 0; idx < num_part_segments; idx++)
289 mesher_input.segmentlist[2 * idx + 2 * num_segments] = int(segments[idx].first + num_points);
290 mesher_input.segmentlist[2 * idx + 1 + 2*num_segments] = int(segments[idx].second + num_points);
296 EXCEPTION(
"Either a part, vtk surface or both are required for 2d meshing");
301 mesher_input.holelist = (
double *) malloc(
mHoles.size() * 2 *
sizeof(double));
302 mesher_input.numberofholes = int(
mHoles.size());
303 for (
unsigned idx = 0; idx <
mHoles.size(); idx++)
306 for (
unsigned jdx = 0; jdx < 2; jdx++)
308 mesher_input.holelist[2 * idx + jdx] = hole_location[jdx];
315 mesher_input.regionlist = (
double *) malloc(
mRegions.size() * 4 *
sizeof(double));
316 mesher_input.numberofregions = int(
mRegions.size());
317 for (
unsigned idx = 0; idx <
mRegions.size(); idx++)
320 for (
unsigned jdx = 0; jdx < 2; jdx++)
322 mesher_input.regionlist[4 * idx + jdx] = region_location[jdx];
324 mesher_input.regionlist[4 * idx + 2] = 1.0;
329 std::string mesher_command =
"pqQze";
333 mesher_command +=
"a" + boost::lexical_cast<std::string>(mesh_size);
339 mesher_command +=
"A";
343 triangulate((
char*) mesher_command.c_str(), &mesher_input, &mesher_output, NULL);
345 mpMesh->ImportFromMesher(mesher_output, mesher_output.numberoftriangles, mesher_output.trianglelist,
346 mesher_output.numberofedges, mesher_output.edgelist, mesher_output.edgemarkerlist);
351 unsigned num_elements = mesher_output.numberoftriangles;
352 for(
unsigned idx=0; idx< num_elements; idx++)
354 mAttributes.push_back(
double(mesher_output.triangleattributelist[idx]));
364 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
367 std::vector<DimensionalChastePoint<SPACE_DIM> > vertex_locations =
mpDomain->GetVertexLocations();
368 std::vector<DimensionalChastePoint<SPACE_DIM> > hole_locations =
mpDomain->GetHoleMarkers();
369 unsigned num_vertices = vertex_locations.size();
370 unsigned num_holes = hole_locations.size();
371 std::vector<boost::shared_ptr<Facet<SPACE_DIM> > > facets =
mpDomain->GetFacets();
372 unsigned num_facets = facets.size();
374 class tetgen::tetgenio mesher_input, mesher_output;
376 tetgen::tetgenio::facet *f;
377 tetgen::tetgenio::polygon *p;
378 mesher_input.pointlist =
new double[(num_vertices) * 3];
379 mesher_input.numberofpoints = num_vertices;
381 for (
unsigned idx = 0; idx < num_vertices; idx++)
383 c_vector<double, SPACE_DIM> vertex_location = vertex_locations[idx].GetLocation(
mReferenceLength);
384 for (
unsigned jdx = 0; jdx < 3; jdx++)
386 mesher_input.pointlist[3 * idx + jdx] = vertex_location[jdx];
391 mesher_input.holelist =
new double[(num_holes) * 3];
392 mesher_input.numberofholes = num_holes;
393 for (
unsigned idx = 0; idx < num_holes; idx++)
395 c_vector<double, SPACE_DIM> hole_location = hole_locations[idx].GetLocation(
mReferenceLength);
396 for (
unsigned jdx = 0; jdx < 3; jdx++)
398 mesher_input.holelist[3 * idx + jdx] = hole_location[jdx];
402 mesher_input.numberoffacets = num_facets;
403 mesher_input.facetlist =
new tetgen::tetgenio::facet[num_facets];
404 mesher_input.facetmarkerlist =
new int[num_facets];
405 for (
unsigned idx = 0; idx < num_facets; idx++)
407 mesher_input.facetmarkerlist[idx] = 0;
408 f = &mesher_input.facetlist[idx];
409 std::vector<boost::shared_ptr<Polygon<SPACE_DIM> > > polygons = facets[idx]->GetPolygons();
410 f->numberofpolygons = polygons.size();
411 f->polygonlist =
new tetgen::tetgenio::polygon[f->numberofpolygons];
412 f->numberofholes = 0;
414 for (
unsigned jdx = 0; jdx < polygons.size(); jdx++)
416 p = &f->polygonlist[jdx];
417 p->numberofvertices = polygons[jdx]->GetVertices().size();
418 p->vertexlist =
new int[p->numberofvertices];
419 for (
unsigned kdx = 0; kdx < polygons[jdx]->GetVertices().size(); kdx++)
421 p->vertexlist[kdx] = int(polygons[jdx]->GetVertices()[kdx]->GetIndex());
428 unsigned num_holes =
mHoles.size();
429 mesher_input.holelist =
new double[(num_holes) * 3];
430 mesher_input.numberofholes = num_holes;
431 for (
unsigned idx = 0; idx < num_holes; idx++)
434 for (
unsigned jdx = 0; jdx < 3; jdx++)
436 mesher_input.holelist[3 * idx + jdx] = hole_location[jdx];
441 std::string mesher_command =
"pqQz";
445 mesher_command +=
"a" + boost::lexical_cast<std::string>(mesh_size);
449 tetgen::tetrahedralize((
char*) mesher_command.c_str(), &mesher_input, &mesher_output);
450 mpMesh->ImportDiscreteContinuumMeshFromTetgen(mesher_output, mesher_output.numberoftetrahedra, mesher_output.tetrahedronlist,
451 mesher_output.numberoftrifaces, mesher_output.trifacelist, NULL);
455 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
458 class tetgen::tetgenio mesher_input, mesher_output;
462 mesher_input.load_stl(writable);
466 unsigned num_holes =
mHoles.size();
467 mesher_input.holelist =
new double[(num_holes) * 3];
468 mesher_input.numberofholes = num_holes;
469 for (
unsigned idx = 0; idx < num_holes; idx++)
472 for (
unsigned jdx = 0; jdx < 3; jdx++)
474 mesher_input.holelist[3 * idx + jdx] = hole_location[jdx];
479 std::string mesher_command =
"pqQz";
483 mesher_command +=
"a" + boost::lexical_cast<std::string>(mesh_size);
487 tetgen::tetrahedralize((
char*) mesher_command.c_str(), &mesher_input, &mesher_output);
489 mpMesh->ImportDiscreteContinuumMeshFromTetgen(mesher_output, mesher_output.numberoftetrahedra, mesher_output.tetrahedronlist,
490 mesher_output.numberoftrifaces, mesher_output.trifacelist, NULL);
495 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
498 mesherIo.numberofpoints = 0;
499 mesherIo.pointlist = NULL;
500 mesherIo.numberofpointattributes = 0;
501 mesherIo.pointattributelist = (
double *) NULL;
502 mesherIo.pointmarkerlist = (
int *) NULL;
503 mesherIo.segmentlist = NULL;
504 mesherIo.segmentmarkerlist = (
int *) NULL;
505 mesherIo.numberofsegments = 0;
506 mesherIo.numberofholes = 0;
507 mesherIo.numberofregions = 0;
508 mesherIo.trianglelist = (
int *) NULL;
509 mesherIo.triangleattributelist = (
double *) NULL;
510 mesherIo.numberoftriangleattributes = 0;
511 mesherIo.edgelist = (
int *) NULL;
512 mesherIo.edgemarkerlist = (
int *) NULL;
515 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
518 if (mesherIo.numberofpoints != 0)
520 mesherIo.numberofpoints = 0;
521 free(mesherIo.pointlist);
524 if (mesherIo.numberofsegments != 0)
526 mesherIo.numberofsegments = 0;
527 free(mesherIo.segmentlist);
531 free(mesherIo.pointattributelist);
532 free(mesherIo.pointmarkerlist);
533 free(mesherIo.segmentmarkerlist);
534 free(mesherIo.trianglelist);
535 free(mesherIo.triangleattributelist);
536 free(mesherIo.edgelist);
537 free(mesherIo.edgemarkerlist);
std::vector< DimensionalChastePoint< SPACE_DIM > > mHoles
~DiscreteContinuumMeshGenerator()
void SetDomain(boost::shared_ptr< Part< SPACE_DIM > > pDomain)
std::vector< unsigned > mAttributes
units::quantity< unit::length > mReferenceLength
void SetRegionMarkers(std::vector< DimensionalChastePoint< SPACE_DIM > > regionMarkers)
void SetMaxElementArea(units::quantity< unit::volume > maxElementArea)
boost::shared_ptr< DiscreteContinuumMesh< ELEMENT_DIM, SPACE_DIM > > mpMesh
static boost::shared_ptr< DiscreteContinuumMeshGenerator< ELEMENT_DIM, SPACE_DIM > > Create()
void InitialiseTriangulateIo(triangulateio &mesherIo)
void FreeTriangulateIo(triangulateio &mesherIo)
boost::shared_ptr< DiscreteContinuumMesh< ELEMENT_DIM, SPACE_DIM > > GetMesh()
void SetHoles(std::vector< DimensionalChastePoint< SPACE_DIM > > holes)
boost::shared_ptr< Part< SPACE_DIM > > mpDomain
vtkSmartPointer< vtkPolyData > mpVtkDomain
std::vector< DimensionalChastePoint< SPACE_DIM > > mRegions
units::quantity< unit::volume > mMaxElementArea
DiscreteContinuumMeshGenerator()