Go to the source code of this file.
◆ add_cube_convex_hull_to_mesh()
      
        
          | void add_cube_convex_hull_to_mesh  | 
          ( | 
          MeshBase &  | 
          mesh,  | 
        
        
           | 
           | 
          Point  | 
          lower_limit,  | 
        
        
           | 
           | 
          Point  | 
          upper_limit  | 
        
        
           | 
          ) | 
           |  | 
        
      
 
Definition at line 249 of file miscellaneous_ex6.C.
  253 #ifdef LIBMESH_HAVE_TETGEN 
  260                                     lower_limit(0), upper_limit(0),
 
  261                                     lower_limit(1), upper_limit(1),
 
  262                                     lower_limit(2), upper_limit(2),
 
  269   t.pointset_convexhull();
 
  275   std::map<unsigned, unsigned> node_id_map;
 
  276   typedef std::map<unsigned, unsigned>::iterator iterator;
 
  278   for (
auto & elem : cube_mesh.element_ptr_range())
 
  279     for (
auto s : elem->side_index_range())
 
  280       if (elem->neighbor_ptr(s) == 
nullptr)
 
  283           std::unique_ptr<Elem> side = elem->side_ptr(s);
 
  285           for (
auto n : side->node_index_range())
 
  286             node_id_map.insert(std::make_pair(side->node_id(n), 0));
 
  291   for (iterator it=node_id_map.begin(); it != node_id_map.end(); ++it)
 
  294       unsigned id = (*it).first;
 
  297       Node & old_node = cube_mesh.node_ref(
id);
 
  303       (*it).second = new_node->
id();
 
  309   for (
auto & old_elem : cube_mesh.element_ptr_range())
 
  310     if (old_elem->type() == 
TRI3)
 
  316         for (
auto i : old_elem->node_index_range())
 
  319             iterator it = node_id_map.find(old_elem->node_id(i));
 
  322             if (it == node_id_map.end())
 
  323               libmesh_error_msg(
"Node id " << old_elem->node_id(i) << 
" not found in map!");
 
  326             unsigned new_node_id = (*it).second;
 
  335 #endif // LIBMESH_HAVE_TETGEN 
 
References libMesh::MeshBase::add_elem(), libMesh::MeshBase::add_point(), libMesh::MeshTools::Generation::build_cube(), libMesh::ParallelObject::comm(), libMesh::HEX8, libMesh::DofObject::id(), libMesh::libmesh_ignore(), mesh, libMesh::MeshTools::n_elem(), libMesh::MeshBase::node_ptr(), libMesh::TetGenMeshInterface::pointset_convexhull(), libMesh::Elem::set_node(), and libMesh::TRI3.
Referenced by tetrahedralize_domain().
 
 
◆ main()
      
        
          | int main  | 
          ( | 
          int  | 
          argc,  | 
        
        
           | 
           | 
          char **  | 
          argv  | 
        
        
           | 
          ) | 
           |  | 
        
      
 
 
◆ tetrahedralize_domain()
      
        
          | void tetrahedralize_domain  | 
          ( | 
          const Parallel::Communicator &  | 
          comm | ) | 
           | 
        
      
 
Definition at line 168 of file miscellaneous_ex6.C.
  170 #ifdef LIBMESH_HAVE_TETGEN 
  182   Point hole_lower_limit(0.2, 0.2, 0.4);
 
  183   Point hole_upper_limit(0.8, 0.8, 0.6);
 
  197   std::vector<Point> hole(1);
 
  198   hole[0] = 
Point(0.5*(hole_lower_limit + hole_upper_limit));
 
  203   double quality_constraint = 2.0;
 
  208   double volume_constraint = 0.001;
 
  220   t.set_switches(
"VCC");
 
  222   t.triangulate_conformingDelaunayMesh_carvehole(hole,
 
  234 #endif // LIBMESH_HAVE_TETGEN 
 
References add_cube_convex_hull_to_mesh(), libMesh::MeshBase::find_neighbors(), libMesh::libmesh_ignore(), mesh, libMesh::MeshBase::prepare_for_use(), libMesh::TetGenMeshInterface::set_switches(), libMesh::TetGenMeshInterface::triangulate_conformingDelaunayMesh_carvehole(), and libMesh::MeshBase::write().
Referenced by main().
 
 
◆ triangulate_domain()
      
        
          | void triangulate_domain  | 
          ( | 
          const Parallel::Communicator &  | 
          comm | ) | 
           | 
        
      
 
Definition at line 79 of file miscellaneous_ex6.C.
   81 #ifdef LIBMESH_HAVE_TRIANGLE 
  104   t.desired_area() = .01;
 
  112   t.triangulation_type() = TriangleInterface::PSLG;
 
  116   t.smooth_after_generating() = 
true;
 
  121   PolygonHole hole_1(
Point(-0.5,  0.5), 
 
  126   PolygonHole hole_2(
Point(0.5, 0.5),   
 
  131   Point ellipse_center(-0.5,  -0.5);
 
  132   const unsigned int n_ellipse_points=100;
 
  133   std::vector<Point> ellipse_points(n_ellipse_points);
 
  135     dtheta = 2*
libMesh::pi / static_cast<Real>(n_ellipse_points),
 
  139   for (
unsigned int i=0; i<n_ellipse_points; ++i)
 
  140     ellipse_points[i]= 
Point(ellipse_center(0)+a*cos(i*dtheta),
 
  141                              ellipse_center(1)+b*sin(i*dtheta));
 
  143   ArbitraryHole hole_3(ellipse_center, ellipse_points);
 
  146   std::vector<Hole *> holes;
 
  147   holes.push_back(&hole_1);
 
  148   holes.push_back(&hole_2);
 
  149   holes.push_back(&hole_3);
 
  152   t.attach_hole_list(&holes);
 
  163 #endif // LIBMESH_HAVE_TRIANGLE 
 
References libMesh::MeshBase::add_point(), libMesh::TriangleInterface::attach_hole_list(), libMesh::TriangleInterface::desired_area(), libMesh::libmesh_ignore(), mesh, libMesh::pi, libMesh::TriangleInterface::PSLG, libMesh::Real, libMesh::TriangleInterface::smooth_after_generating(), libMesh::TriangleInterface::triangulate(), libMesh::TriangleInterface::triangulation_type(), and libMesh::MeshBase::write().
Referenced by main().
 
 
 
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
 
A concrete instantiation of the Hole class that describes polygonal (triangular, square,...
 
virtual void write(const std::string &name)=0
 
const Parallel::Communicator & comm() const
 
virtual const Node * node_ptr(const dof_id_type i) const =0
 
void tetrahedralize_domain(const Parallel::Communicator &comm)
 
An abstract class for defining a 2-dimensional hole.
 
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
 
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
 
void triangulate_domain(const Parallel::Communicator &comm)
 
void libmesh_ignore(const Args &...)
 
A Point defines a location in LIBMESH_DIM dimensional Real space.
 
A Node is like a Point, but with more information.
 
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
 
virtual Node *& set_node(const unsigned int i)
 
Class TetGenMeshInterface provides an interface for tetrahedralization of meshes using the TetGen lib...
 
void add_cube_convex_hull_to_mesh(MeshBase &mesh, Point lower_limit, Point upper_limit)
 
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
 
The Tri3 is an element in 2D composed of 3 nodes.
 
This is the base class from which all geometric element types are derived.
 
Another concrete instantiation of the hole, this one should be sufficiently general for most non-poly...
 
A C++ interface between LibMesh and the Triangle library written by J.R.
 
virtual Node * add_point(const Point &p, const dof_id_type id=DofObject::invalid_id, const processor_id_type proc_id=DofObject::invalid_processor_id)=0
Add a new Node at Point p to the end of the vertex array, with processor_id procid.
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
void prepare_for_use(const bool skip_renumber_nodes_and_elements=false, const bool skip_find_neighbors=false)
Prepare a newly ecreated (or read) mesh for use.
 
virtual void find_neighbors(const bool reset_remote_elements=false, const bool reset_current_list=true)=0
Locate element face (edge in 2D) neighbors.