$extrastylesheet
00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 #include "libmesh/libmesh_config.h" 00020 00021 #ifdef LIBMESH_HAVE_TRIANGLE 00022 00023 // C/C++ includes 00024 #include <sstream> 00025 00026 #include "libmesh/mesh_triangle_interface.h" 00027 #include "libmesh/unstructured_mesh.h" 00028 #include "libmesh/face_tri3.h" 00029 #include "libmesh/face_tri6.h" 00030 #include "libmesh/mesh_generation.h" 00031 #include "libmesh/mesh_smoother_laplace.h" 00032 #include "libmesh/boundary_info.h" 00033 #include "libmesh/mesh_triangle_holes.h" 00034 #include "libmesh/mesh_triangle_wrapper.h" 00035 00036 namespace libMesh 00037 { 00038 // 00039 // Function definitions for the TriangleInterface class 00040 // 00041 00042 // Constructor 00043 TriangleInterface::TriangleInterface(UnstructuredMesh& mesh) 00044 : _mesh(mesh), 00045 _holes(NULL), 00046 _elem_type(TRI3), 00047 _desired_area(0.1), 00048 _minimum_angle(20.0), 00049 _extra_flags(""), 00050 _triangulation_type(GENERATE_CONVEX_HULL), 00051 _insert_extra_points(false), 00052 _smooth_after_generating(true), 00053 _serializer(_mesh) 00054 {} 00055 00056 00057 00058 // Primary function responsible for performing the triangulation 00059 void TriangleInterface::triangulate() 00060 { 00061 // Will the triangulation have holes? 00062 const bool have_holes = ((_holes != NULL) && (!_holes->empty())); 00063 00064 // If the initial PSLG is really simple, e.g. an L-shaped domain or 00065 // a square/rectangle, the resulting triangulation may be very 00066 // "structured" looking. Sometimes this is a problem if your 00067 // intention is to work with an "unstructured" looking grid. We can 00068 // attempt to work around this limitation by inserting midpoints 00069 // into the original PSLG. Inserting additional points into a 00070 // set of points meant to be a convex hull usually makes less sense. 00071 00072 // May or may not need to insert new points ... 00073 if ((_triangulation_type==PSLG) && (_insert_extra_points)) 00074 { 00075 // Make a copy of the original points from the Mesh 00076 std::vector<Point> original_points (_mesh.n_nodes()); 00077 00078 MeshBase::node_iterator node_it = _mesh.nodes_begin(); 00079 const MeshBase::node_iterator node_end = _mesh.nodes_end(); 00080 00081 for (unsigned int ctr=0; node_it != node_end; ++node_it) 00082 original_points[ctr++] = **node_it; 00083 00084 // Clear out the mesh 00085 _mesh.clear(); 00086 00087 // Make sure the new Mesh will be 2D 00088 _mesh.set_mesh_dimension(2); 00089 00090 // Insert a new point on each PSLG at some random location 00091 // np=index into new points vector 00092 // n =index into original points vector 00093 for (unsigned int np=0, n=0; np<2*original_points.size(); ++np) 00094 { 00095 // the even entries are the original points 00096 if (np%2==0) 00097 _mesh.add_point(original_points[n++]); 00098 00099 else // the odd entries are the midpoints of the original PSLG segments 00100 _mesh.add_point ((original_points[n] + original_points[n-1])/2); 00101 } 00102 } 00103 00104 // Regardless of whether we added additional points, the set of points to 00105 // triangulate is now sitting in the mesh. 00106 00107 // If the holes vector is non-NULL (and non-empty) we need to determine 00108 // the number of additional points which the holes will add to the 00109 // triangulation. 00110 unsigned int n_hole_points = 0; 00111 00112 if (have_holes) 00113 { 00114 for (unsigned int i=0; i<_holes->size(); ++i) 00115 n_hole_points += (*_holes)[i]->n_points(); 00116 } 00117 00118 // Triangle data structure for the mesh 00119 TriangleWrapper::triangulateio initial; 00120 TriangleWrapper::triangulateio final; 00121 00122 // Pseudo-Constructor for the triangle io structs 00123 TriangleWrapper::init(initial); 00124 TriangleWrapper::init(final); 00125 00126 initial.numberofpoints = _mesh.n_nodes() + n_hole_points; 00127 initial.pointlist = static_cast<REAL*>(std::malloc(initial.numberofpoints * 2 * sizeof(REAL))); 00128 00129 if (_triangulation_type==PSLG) 00130 { 00131 // Implicit segment ordering: One segment per point, including hole points 00132 if (this->segments.empty()) 00133 initial.numberofsegments = initial.numberofpoints; 00134 00135 // User-defined segment ordering: One segment per entry in the segments vector 00136 else 00137 initial.numberofsegments = this->segments.size(); 00138 } 00139 00140 else if (_triangulation_type==GENERATE_CONVEX_HULL) 00141 initial.numberofsegments = n_hole_points; // One segment for each hole point 00142 00143 // Debugging 00144 // libMesh::out << "Number of segments set to: " << initial.numberofsegments << std::endl; 00145 00146 // Allocate space for the segments (2 int per segment) 00147 if (initial.numberofsegments > 0) 00148 { 00149 initial.segmentlist = static_cast<int*> (std::malloc(initial.numberofsegments * 2 * sizeof(int))); 00150 } 00151 00152 00153 // Copy all the holes' points and segments into the triangle struct. 00154 00155 // The hole_offset is a constant offset into the points vector which points 00156 // past the end of the last hole point added. 00157 unsigned int hole_offset=0; 00158 00159 if (have_holes) 00160 for (unsigned int i=0; i<_holes->size(); ++i) 00161 { 00162 for (unsigned int ctr=0, h=0; h<(*_holes)[i]->n_points(); ctr+=2, ++h) 00163 { 00164 Point p = (*_holes)[i]->point(h); 00165 00166 const unsigned int index0 = 2*hole_offset+ctr; 00167 const unsigned int index1 = 2*hole_offset+ctr+1; 00168 00169 // Save the x,y locations in the triangle struct. 00170 initial.pointlist[index0] = p(0); 00171 initial.pointlist[index1] = p(1); 00172 00173 // Set the points which define the segments 00174 initial.segmentlist[index0] = hole_offset+h; 00175 initial.segmentlist[index1] = (h==(*_holes)[i]->n_points()-1) ? hole_offset : hole_offset+h+1; // wrap around 00176 } 00177 00178 // Update the hole_offset for the next hole 00179 hole_offset += (*_holes)[i]->n_points(); 00180 } 00181 00182 00183 // Copy all the non-hole points and segments into the triangle struct. 00184 { 00185 MeshBase::node_iterator it = _mesh.nodes_begin(); 00186 const MeshBase::node_iterator end = _mesh.nodes_end(); 00187 00188 for (dof_id_type ctr=0; it != end; ctr+=2, ++it) 00189 { 00190 dof_id_type index = 2*hole_offset + ctr; 00191 00192 // Get pointer to the current node 00193 Node* node = *it; 00194 00195 // Set x,y values in pointlist 00196 initial.pointlist[index] = (*node)(0); 00197 initial.pointlist[index+1] = (*node)(1); 00198 00199 // If the user requested a PSLG, the non-hole points are also segments 00200 if (_triangulation_type==PSLG) 00201 { 00202 // Use implicit ordering to define segments 00203 if (this->segments.empty()) 00204 { 00205 dof_id_type n = ctr/2; // ctr is always even 00206 initial.segmentlist[index] = hole_offset+n; 00207 initial.segmentlist[index+1] = (n==_mesh.n_nodes()-1) ? hole_offset : hole_offset+n+1; // wrap around 00208 } 00209 } 00210 } 00211 } 00212 00213 00214 // If the user provided it, use his ordering to define the segments 00215 for (unsigned int ctr=0, s=0; s<this->segments.size(); ctr+=2, ++s) 00216 { 00217 const unsigned int index0 = 2*hole_offset+ctr; 00218 const unsigned int index1 = 2*hole_offset+ctr+1; 00219 00220 initial.segmentlist[index0] = hole_offset + this->segments[s].first; 00221 initial.segmentlist[index1] = hole_offset + this->segments[s].second; 00222 } 00223 00224 00225 00226 // Tell the input struct about the holes 00227 if (have_holes) 00228 { 00229 initial.numberofholes = _holes->size(); 00230 initial.holelist = static_cast<REAL*>(std::malloc(initial.numberofholes * 2 * sizeof(REAL))); 00231 for (unsigned int i=0, ctr=0; i<_holes->size(); ++i, ctr+=2) 00232 { 00233 Point inside_point = (*_holes)[i]->inside(); 00234 initial.holelist[ctr] = inside_point(0); 00235 initial.holelist[ctr+1] = inside_point(1); 00236 } 00237 } 00238 00239 // Set the triangulation flags. 00240 // c ~ enclose convex hull with segments 00241 // z ~ use zero indexing 00242 // B ~ Suppresses boundary markers in the output 00243 // Q ~ run in "quiet" mode 00244 // p ~ Triangulates a Planar Straight Line Graph 00245 // If the `p' switch is used, `segmentlist' must point to a list of 00246 // segments, `numberofsegments' must be properly set, and 00247 // `segmentmarkerlist' must either be set to NULL (in which case all 00248 // markers default to zero), or must point to a list of markers. 00249 // D ~ Conforming Delaunay: use this switch if you want all triangles 00250 // in the mesh to be Delaunay, and not just constrained Delaunay 00251 // q ~ Quality mesh generation with no angles smaller than 20 degrees. 00252 // An alternate minimum angle may be specified after the q 00253 // a ~ Imposes a maximum triangle area constraint. 00254 // -P Suppresses the output .poly file. Saves disk space, but you lose the ability to maintain 00255 // constraining segments on later refinements of the mesh. 00256 // Create the flag strings, depends on element type 00257 std::ostringstream flags; 00258 00259 // Default flags always used 00260 flags << "zBPQ"; 00261 00262 // Flags which are specific to the type of triangulation 00263 switch (_triangulation_type) 00264 { 00265 case GENERATE_CONVEX_HULL: 00266 { 00267 flags << "c"; 00268 break; 00269 } 00270 00271 case PSLG: 00272 { 00273 flags << "p"; 00274 break; 00275 } 00276 00277 case INVALID_TRIANGULATION_TYPE: 00278 libmesh_error_msg("ERROR: INVALID_TRIANGULATION_TYPE selected!"); 00279 00280 default: 00281 libmesh_error_msg("Unrecognized _triangulation_type"); 00282 } 00283 00284 00285 // Flags specific to the type of element 00286 switch (_elem_type) 00287 { 00288 case TRI3: 00289 { 00290 // do nothing. 00291 break; 00292 } 00293 00294 case TRI6: 00295 { 00296 flags << "o2"; 00297 break; 00298 } 00299 00300 default: 00301 libmesh_error_msg("ERROR: Unrecognized triangular element type."); 00302 } 00303 00304 00305 // If we do have holes and the user asked to GENERATE_CONVEX_HULL, 00306 // need to add the p flag so the triangulation respects those segments. 00307 if ((_triangulation_type==GENERATE_CONVEX_HULL) && (have_holes)) 00308 flags << "p"; 00309 00310 // Finally, add the area constraint 00311 if (_desired_area > TOLERANCE) 00312 flags << "a" << std::fixed << _desired_area; 00313 00314 // add minimum angle constraint 00315 if (_minimum_angle > TOLERANCE) 00316 flags << "q" << std::fixed << _minimum_angle; 00317 00318 // add user provided extra flags 00319 if (_extra_flags.size() > 0) 00320 flags << _extra_flags; 00321 00322 // Refine the initial output to conform to the area constraint 00323 TriangleWrapper::triangulate(const_cast<char*>(flags.str().c_str()), 00324 &initial, 00325 &final, 00326 NULL); // voronoi ouput -- not used 00327 00328 00329 // Send the information computed by Triangle to the Mesh. 00330 TriangleWrapper::copy_tri_to_mesh(final, 00331 _mesh, 00332 _elem_type); 00333 00334 // To the naked eye, a few smoothing iterations usually looks better, 00335 // so we do this by default unless the user says not to. 00336 if (this->_smooth_after_generating) 00337 LaplaceMeshSmoother(_mesh).smooth(2); 00338 00339 00340 // Clean up. 00341 TriangleWrapper::destroy(initial, TriangleWrapper::INPUT); 00342 TriangleWrapper::destroy(final, TriangleWrapper::OUTPUT); 00343 00344 } 00345 00346 } // namespace libMesh 00347 00348 00349 00350 00351 00352 00353 00354 #endif // LIBMESH_HAVE_TRIANGLE