Chaste  Build::
DiscreteContinuumMeshGenerator.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 <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>
40 #include <vtkLine.h>
41 #include <vtkPoints.h>
42 #include "Exception.hpp"
43 #include "Warnings.hpp"
44 #include "Facet.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"
53 
54 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
56  mMaxElementArea(0.0 * unit::metres * unit::metres* unit::metres),
57  mpMesh(),
58  mpDomain(),
59  mpVtkDomain(),
60  mStlFilePath(),
61  mHoles(),
62  mRegions(),
63  mAttributes(),
64  mReferenceLength(BaseUnits::Instance()->GetReferenceLengthScale())
65 {
66 
67 }
68 
69 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
70 boost::shared_ptr<DiscreteContinuumMeshGenerator<ELEMENT_DIM, SPACE_DIM> > DiscreteContinuumMeshGenerator<ELEMENT_DIM, SPACE_DIM>::Create()
71 {
73  return pSelf;
74 }
75 
76 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
78 {
79 
80 }
81 
82 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
83 boost::shared_ptr<DiscreteContinuumMesh<ELEMENT_DIM, SPACE_DIM> > DiscreteContinuumMeshGenerator<ELEMENT_DIM, SPACE_DIM>::GetMesh()
84 {
85  if(!mpMesh)
86  {
87  EXCEPTION("No mesh has been generated.");
88  }
89  return mpMesh;
90 }
91 
92 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
94 {
95  mpDomain = pDomain;
96 }
97 
98 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
100 {
101  mpVtkDomain = pDomain;
102 }
103 
104 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
106 {
107  mStlFilePath = rPathToStl;
108 }
109 
110 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
111 void DiscreteContinuumMeshGenerator<ELEMENT_DIM, SPACE_DIM>::SetMaxElementArea(units::quantity<unit::volume> maxElementArea)
112 {
113  mMaxElementArea = maxElementArea;
114 }
115 
116 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
118 {
119  mHoles = holes;
120 }
121 
122 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
124 {
125  mRegions = regionMarkers;
126 }
127 
128 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
130 {
131  // Create a mesh
132  mpMesh = boost::shared_ptr<DiscreteContinuumMesh<ELEMENT_DIM, SPACE_DIM> >(new DiscreteContinuumMesh<ELEMENT_DIM, SPACE_DIM>());
133 
134  // For 2D parts use triangle
135  if (ELEMENT_DIM == 2)
136  {
137  if (SPACE_DIM == 2)
138  {
139  Mesh2d();
140  }
141  else
142  {
143  if(mpDomain)
144  {
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)
147  {
148  Mesh2d();
149  }
150  else
151  {
152  EXCEPTION("2D meshing is only supported for parts with z=0.");
153  }
154  }
155  else if(mpVtkDomain)
156  {
157  Mesh2d();
158  }
159  else
160  {
161  EXCEPTION("2d meshing is not supported for STL files.");
162  }
163  }
164  }
165  // For 3d use tetgen
166  else
167  {
168  if(mpDomain)
169  {
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)
172  {
173  EXCEPTION("The part is two-dimensional, use the 2D meshing functionality.");
174  }
175  else
176  {
177  Mesh3d();
178  }
179  }
180  else if(mpVtkDomain)
181  {
182  Mesh3d();
183  }
184  else
185  {
186  MeshStl3d();
187  }
188  }
189 }
190 
191 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
193 {
194  struct triangulateio mesher_input, mesher_output;
195  this->InitialiseTriangulateIo(mesher_input);
196  this->InitialiseTriangulateIo(mesher_output);
197 
198  // Cases: have just domain, have just vtk domain, have vtk domain and domain
199  if(mpDomain and !mpVtkDomain)
200  {
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++)
206  {
207  c_vector<double, SPACE_DIM> vertex_location = vertex_locations[idx].GetLocation(mReferenceLength);
208  for (unsigned jdx = 0; jdx < 2; jdx++)
209  {
210  mesher_input.pointlist[2 * idx + jdx] = vertex_location[jdx];
211  }
212  }
213 
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++)
219  {
220  mesher_input.segmentlist[2 * idx] = int(segments[idx].first);
221  mesher_input.segmentlist[2 * idx + 1] = int(segments[idx].second);
222  }
223  }
224  else if(!mpDomain and mpVtkDomain)
225  {
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++)
230  {
231  for (unsigned jdx = 0; jdx < 2; jdx++)
232  {
233  mesher_input.pointlist[2 * idx + jdx] = mpVtkDomain->GetPoints()->GetPoint(idx)[jdx];
234  }
235  }
236 
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++)
243  {
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));
247  }
248  }
249  else if(mpDomain and mpVtkDomain)
250  {
251  unsigned num_points = mpVtkDomain->GetNumberOfPoints();
252  std::vector<DimensionalChastePoint<SPACE_DIM> > vertex_locations = mpDomain->GetVertexLocations();
253  unsigned num_vertices = vertex_locations.size();
254 
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++)
258  {
259  for (unsigned jdx = 0; jdx < 2; jdx++)
260  {
261  mesher_input.pointlist[2 * idx + jdx] = mpVtkDomain->GetPoints()->GetPoint(idx)[jdx];
262  }
263  }
264  for (unsigned idx = 0; idx < num_vertices; idx++)
265  {
266  c_vector<double, SPACE_DIM> vertex_location = vertex_locations[idx].GetLocation(mReferenceLength);
267  for (unsigned jdx = 0; jdx < 2; jdx++)
268  {
269  mesher_input.pointlist[2 * idx + jdx + 2*num_points] = vertex_location[jdx];
270  }
271  }
272 
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();
277 
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++)
282  {
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));
286  }
287  for (unsigned idx = 0; idx < num_part_segments; idx++)
288  {
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);
291  }
292 
293  }
294  else
295  {
296  EXCEPTION("Either a part, vtk surface or both are required for 2d meshing");
297  }
298 
299  if (mHoles.size() > 0)
300  {
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++)
304  {
305  c_vector<double, SPACE_DIM> hole_location = mHoles[idx].GetLocation(mReferenceLength);
306  for (unsigned jdx = 0; jdx < 2; jdx++)
307  {
308  mesher_input.holelist[2 * idx + jdx] = hole_location[jdx];
309  }
310  }
311  }
312 
313  if (mRegions.size() > 0)
314  {
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++)
318  {
319  c_vector<double, SPACE_DIM> region_location = mRegions[idx].GetLocation(mReferenceLength);
320  for (unsigned jdx = 0; jdx < 2; jdx++)
321  {
322  mesher_input.regionlist[4 * idx + jdx] = region_location[jdx];
323  }
324  mesher_input.regionlist[4 * idx + 2] = 1.0;
325  mesher_input.regionlist[4 * idx + 3] = mMaxElementArea/units::pow<3>(mReferenceLength);
326  }
327  }
328 
329  std::string mesher_command = "pqQze";
330  if (mMaxElementArea > 0.0*unit::metres*unit::metres*unit::metres)
331  {
332  double mesh_size = mMaxElementArea/units::pow<3>(mReferenceLength);
333  mesher_command += "a" + boost::lexical_cast<std::string>(mesh_size);
334  }
335  if(mRegions.size()>0)
336  {
337  if(mRegions.size()>0)
338  {
339  mesher_command += "A";
340  }
341  }
342 
343  triangulate((char*) mesher_command.c_str(), &mesher_input, &mesher_output, NULL);
344 
345  mpMesh->ImportFromMesher(mesher_output, mesher_output.numberoftriangles, mesher_output.trianglelist,
346  mesher_output.numberofedges, mesher_output.edgelist, mesher_output.edgemarkerlist);
347 
348  mAttributes.clear();
349  if(mRegions.size()>0)
350  {
351  unsigned num_elements = mesher_output.numberoftriangles;
352  for(unsigned idx=0; idx< num_elements; idx++)
353  {
354  mAttributes.push_back(double(mesher_output.triangleattributelist[idx]));
355  }
356  mpMesh->SetAttributes(mAttributes);
357  }
358 
359  //Tidy up triangle
360  this->FreeTriangulateIo(mesher_input);
361  this->FreeTriangulateIo(mesher_output);
362 }
363 
364 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
366 {
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();
373 
374  class tetgen::tetgenio mesher_input, mesher_output;
375 
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;
380 
381  for (unsigned idx = 0; idx < num_vertices; idx++)
382  {
383  c_vector<double, SPACE_DIM> vertex_location = vertex_locations[idx].GetLocation(mReferenceLength);
384  for (unsigned jdx = 0; jdx < 3; jdx++)
385  {
386  mesher_input.pointlist[3 * idx + jdx] = vertex_location[jdx];
387  }
388  }
389 
390  // Add the holes
391  mesher_input.holelist = new double[(num_holes) * 3];
392  mesher_input.numberofholes = num_holes;
393  for (unsigned idx = 0; idx < num_holes; idx++)
394  {
395  c_vector<double, SPACE_DIM> hole_location = hole_locations[idx].GetLocation(mReferenceLength);
396  for (unsigned jdx = 0; jdx < 3; jdx++)
397  {
398  mesher_input.holelist[3 * idx + jdx] = hole_location[jdx];
399  }
400  }
401 
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++)
406  {
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;
413  f->holelist = NULL;
414  for (unsigned jdx = 0; jdx < polygons.size(); jdx++)
415  {
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++)
420  {
421  p->vertexlist[kdx] = int(polygons[jdx]->GetVertices()[kdx]->GetIndex());
422  }
423  }
424  }
425 
426  if(mHoles.size() > 0)
427  {
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++)
432  {
433  c_vector<double, SPACE_DIM> hole_location = mHoles[idx].GetLocation(mReferenceLength);
434  for (unsigned jdx = 0; jdx < 3; jdx++)
435  {
436  mesher_input.holelist[3 * idx + jdx] = hole_location[jdx];
437  }
438  }
439  }
440 
441  std::string mesher_command = "pqQz";
442  if (mMaxElementArea > 0.0*unit::metres*unit::metres*unit::metres)
443  {
444  double mesh_size = mMaxElementArea/units::pow<3>(mReferenceLength);
445  mesher_command += "a" + boost::lexical_cast<std::string>(mesh_size);
446  }
447 
448  // Library call
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);
452 
453 }
454 
455 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
457 {
458  class tetgen::tetgenio mesher_input, mesher_output;
459  char * writable = new char[mStlFilePath.size() + 1];
460  std::copy(mStlFilePath.begin(), mStlFilePath.end(), writable);
461  writable[mStlFilePath.size()] = '\0';
462  mesher_input.load_stl(writable);
463 
464  if(mHoles.size() > 0)
465  {
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++)
470  {
471  c_vector<double, SPACE_DIM> hole_location = mHoles[idx].GetLocation(mReferenceLength);
472  for (unsigned jdx = 0; jdx < 3; jdx++)
473  {
474  mesher_input.holelist[3 * idx + jdx] = hole_location[jdx];
475  }
476  }
477  }
478 
479  std::string mesher_command = "pqQz";
480  if (mMaxElementArea > 0.0*unit::metres*unit::metres*unit::metres)
481  {
482  double mesh_size = mMaxElementArea/units::pow<3>(mReferenceLength);
483  mesher_command += "a" + boost::lexical_cast<std::string>(mesh_size);
484  }
485 
486  // Library call
487  tetgen::tetrahedralize((char*) mesher_command.c_str(), &mesher_input, &mesher_output);
488 
489  mpMesh->ImportDiscreteContinuumMeshFromTetgen(mesher_output, mesher_output.numberoftetrahedra, mesher_output.tetrahedronlist,
490  mesher_output.numberoftrifaces, mesher_output.trifacelist, NULL);
491 
492  delete[] writable;
493 }
494 
495 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
497 {
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;
513 }
514 
515 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
517 {
518  if (mesherIo.numberofpoints != 0)
519  {
520  mesherIo.numberofpoints = 0;
521  free(mesherIo.pointlist);
522  }
523 
524  if (mesherIo.numberofsegments != 0)
525  {
526  mesherIo.numberofsegments = 0;
527  free(mesherIo.segmentlist);
528  }
529 
530  // These (and the above) should actually be safe since we explicity set to NULL above
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);
538 }
539 
540 // Explicit instantiation
541 template class DiscreteContinuumMeshGenerator<2> ;
542 template class DiscreteContinuumMeshGenerator<3> ;
std::vector< DimensionalChastePoint< SPACE_DIM > > mHoles
void SetDomain(boost::shared_ptr< Part< SPACE_DIM > > pDomain)
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
Definition: Part.hpp:60