quad_from_triangle_mesh.template.h
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2023 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 // Driver for 2D moving block
27 #ifndef OOMPH_QUAD_FROM_TRIANGLE_MESH_TEMPLATE_HEADER
28 #define OOMPH_QUAD_FROM_TRIANGLE_MESH_TEMPLATE_HEADER
29 
30 
31 #ifdef OOMPH_HAS_MPI
32 // MPI header
33 #include "mpi.h"
34 #endif
35 
36 // Standards
37 #include <float.h>
38 #include <iostream>
39 #include <fstream>
40 #include <string.h>
41 #include <iomanip>
42 
43 #ifdef OOMPH_HAS_FPUCONTROLH
44 #include <fpu_control.h>
45 #endif
46 
47 // The mesh
48 #include "../generic/problem.h"
49 #include "../generic/quad_mesh.h"
50 #include "triangle_mesh.template.h"
51 #include "../generic/triangle_scaffold_mesh.h"
52 #include "../generic/unstructured_two_d_mesh_geometry_base.h"
53 #include "../generic/refineable_quad_mesh.h"
54 #include "../generic/Qelements.h"
55 
56 
57 /// /////////////////////////////////////////////////////////////////
58 /// /////////////////////////////////////////////////////////////////
59 /// /////////////////////////////////////////////////////////////////
60 
61 
62 namespace oomph
63 {
64  //============start_of_quad_triangle_class==============================
65  /// Quad mesh built on top of triangle scaffold mesh coming
66  /// from the triangle mesh generator Triangle.
67  /// http://www.cs.cmu.edu/~quake/triangle.html
68  //======================================================================
69  template<class ELEMENT>
71  public virtual QuadMeshBase
72  {
73  public:
74  /// Empty constructor
76  {
77 #ifdef OOMPH_HAS_TRIANGLE_LIB
78  // By default allow the automatic creation of vertices along the
79  // boundaries by Triangle
81 #endif
82 
83  // Mesh can only be built with 2D Qelements.
84  MeshChecker::assert_geometric_element<QElementGeometricBase, ELEMENT>(2);
85  }
86 
87  /// Constructor with the input files
89  const std::string& node_file_name,
90  const std::string& element_file_name,
91  const std::string& poly_file_name,
92  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper,
93  const bool& use_attributes = false,
94  const bool& allow_automatic_creation_of_vertices_on_boundaries = true)
95  {
96  // Mesh can only be built with 2D Qelements.
97  MeshChecker::assert_geometric_element<QElementGeometricBase, ELEMENT>(2);
98 
99  // Initialize the value for allowing creation of points on boundaries
101  allow_automatic_creation_of_vertices_on_boundaries;
102 
103  // Store Timestepper used to build elements
104  this->Time_stepper_pt = time_stepper_pt;
105 
106  // Store the attributes
107  this->Use_attributes = use_attributes;
108 
109  // Build scaffold
110  TriangleScaffoldMesh* tmp_mesh_pt = new TriangleScaffoldMesh(
111  node_file_name, element_file_name, poly_file_name);
112 
113  // Convert mesh from scaffold to actual mesh
114  this->build_from_scaffold(tmp_mesh_pt, time_stepper_pt, use_attributes);
115 
116  // Kill the scaffold
117  delete tmp_mesh_pt;
118  tmp_mesh_pt = 0;
119 
120  // Setup boundary coordinates for boundaries
121  unsigned nbound = nboundary();
122  for (unsigned ibound = 0; ibound < nbound; ibound++)
123  {
124  this->template setup_boundary_coordinates<ELEMENT>(ibound);
125  }
126  }
127 
128 #ifdef OOMPH_HAS_TRIANGLE_LIB
129 
130  /// Build mesh, based on the specifications on
131  /// TriangleMeshParameters. All the actual work is done
132  /// in UnstructuredTwoDMeshGeometryBase
134  TriangleMeshParameters& triangle_mesh_parameters,
135  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper)
136  {
137  // Mesh can only be built with 2D Qelements.
138  MeshChecker::assert_geometric_element<QElementGeometricBase, ELEMENT>(2);
139 
140  // Initialize the value for allowing creation of points on boundaries
142  triangle_mesh_parameters
144 
145  // Store Timestepper used to build elements
146  this->Time_stepper_pt = time_stepper_pt;
147 
148  // ********************************************************************
149  // First part - Get polylines representations
150  // ********************************************************************
151 
152  // Create the polyline representation of all the boundaries and
153  // then create the mesh by calling to "generic_constructor()"
154 
155  // Initialise highest boundary id
156  unsigned max_boundary_id = 0;
157 
158  // *****************************************************************
159  // Part 1.1 - Outer boundary
160  // *****************************************************************
161  // Get the representation of the outer boundaries from the
162  // TriangleMeshParameters object
163  Vector<TriangleMeshClosedCurve*> outer_boundary_pt =
164  triangle_mesh_parameters.outer_boundary_pt();
165 
166 #ifdef PARANOID
167 
168  // Verify that the outer_boundary_object_pt has been set
169  if (outer_boundary_pt.size() == 0)
170  {
171  std::stringstream error_message;
172  error_message
173  << "There are no outer boundaries defined.\n"
174  << "Verify that you have specified the outer boundaries in the\n"
175  << "Triangle_mesh_parameter object\n\n";
176  throw OomphLibError(error_message.str(),
177  OOMPH_CURRENT_FUNCTION,
178  OOMPH_EXCEPTION_LOCATION);
179  } // if (outer_boundary_pt!=0)
180 
181 #endif
182 
183  // Find the number of outer closed curves
184  unsigned n_outer_boundaries = outer_boundary_pt.size();
185 
186  // Create the storage for the polygons that define the outer boundaries
187  Vector<TriangleMeshPolygon*> outer_boundary_polygon_pt(
188  n_outer_boundaries);
189 
190  // Loop over the number of outer boundaries
191  for (unsigned i = 0; i < n_outer_boundaries; ++i)
192  {
193  // Get the polygon representation and compute the max boundary_id on
194  // each outer polygon. Does nothing (i.e. just returns a pointer to
195  // the outer boundary that was input) if the outer boundary is
196  // already a polygon
197  outer_boundary_polygon_pt[i] = this->closed_curve_to_polygon_helper(
198  outer_boundary_pt[i], max_boundary_id);
199  }
200 
201  // *****************************************************************
202  // Part 1.2 - Internal closed boundaries (possible holes)
203  // *****************************************************************
204  // Get the representation of the internal closed boundaries from the
205  // TriangleMeshParameters object
206  Vector<TriangleMeshClosedCurve*> internal_closed_curve_pt =
207  triangle_mesh_parameters.internal_closed_curve_pt();
208 
209  // Find the number of internal closed curves
210  unsigned n_internal_closed_curves = internal_closed_curve_pt.size();
211 
212  // Create the storage for the polygons that define the internal closed
213  // boundaries (again nothing happens (as above) if an internal closed
214  // curve is already a polygon)
215  Vector<TriangleMeshPolygon*> internal_polygon_pt(
216  n_internal_closed_curves);
217 
218  // Loop over the number of internal closed curves
219  for (unsigned i = 0; i < n_internal_closed_curves; ++i)
220  {
221  // Get the polygon representation and compute the max boundary_id on
222  // each internal polygon
223  internal_polygon_pt[i] = this->closed_curve_to_polygon_helper(
224  internal_closed_curve_pt[i], max_boundary_id);
225  }
226 
227  // *****************************************************************
228  // Part 1.3 - Internal open boundaries
229  // *****************************************************************
230  // Get the representation of open boundaries from the
231  // TriangleMeshParameteres object
232  Vector<TriangleMeshOpenCurve*> internal_open_curve_pt =
233  triangle_mesh_parameters.internal_open_curves_pt();
234 
235  // Find the number of internal open curves
236  unsigned n_internal_open_curves = internal_open_curve_pt.size();
237 
238  // Create the storage for the polylines that define the open boundaries
239  Vector<TriangleMeshOpenCurve*> internal_open_curve_poly_pt(
240  n_internal_open_curves);
241 
242  // Loop over the number of internal open curves
243  for (unsigned i = 0; i < n_internal_open_curves; i++)
244  {
245  // Get the open polyline representation and compute the max boundary_id
246  // on each open polyline (again, nothing happens if there are curve
247  // sections on the current internal open curve)
248  internal_open_curve_poly_pt[i] =
250  internal_open_curve_pt[i], max_boundary_id);
251  }
252 
253  // ********************************************************************
254  // Second part - Get associated geom objects and coordinate limits
255  // ********************************************************************
256 
257  // ***************************************************************
258  // Part 2.1 Outer boundary
259  // ***************************************************************
260  for (unsigned i = 0; i < n_outer_boundaries; i++)
261  {
263  outer_boundary_pt[i]);
264  }
265 
266  // ***************************************************************
267  // Part 2.2 - Internal closed boundaries (possible holes)
268  // ***************************************************************
269  for (unsigned i = 0; i < n_internal_closed_curves; i++)
270  {
272  internal_closed_curve_pt[i]);
273  }
274 
275  // ********************************************************************
276  // Part 2.3 - Internal open boundaries
277  // ********************************************************************
278  for (unsigned i = 0; i < n_internal_open_curves; i++)
279  {
281  internal_open_curve_pt[i]);
282  }
283 
284  // ********************************************************************
285  // Third part - Creates the TriangulateIO object by calling the
286  // "generic_constructor()" function
287  // ********************************************************************
288  // Get all the other parameters from the TriangleMeshParameters object
289  // The maximum element area
290  const double element_area = triangle_mesh_parameters.element_area();
291 
292  // The holes coordinates
293  Vector<Vector<double>> extra_holes_coordinates =
294  triangle_mesh_parameters.extra_holes_coordinates();
295 
296  // The regions coordinates
297  std::map<unsigned, Vector<double>> regions =
298  triangle_mesh_parameters.regions_coordinates();
299 
300  // If we use regions then we use attributes
301  const bool use_attributes = triangle_mesh_parameters.is_use_attributes();
302 
303  const bool refine_boundary =
304  triangle_mesh_parameters.is_boundary_refinement_allowed();
305 
306  const bool refine_internal_boundary =
307  triangle_mesh_parameters.is_internal_boundary_refinement_allowed();
308 
309  if (!refine_internal_boundary && refine_boundary)
310  {
311  std::ostringstream error_stream;
312  error_stream
313  << "You have specified that Triangle may refine the outer boundary, "
314  "but\n"
315  << "not internal boundaries. Triangle does not support this "
316  "combination.\n"
317  << "If you do not want Triangle to refine internal boundaries, it "
318  "can't\n"
319  << "refine outer boundaries either!\n"
320  << "Please either disable all boundary refinement\n"
321  << "(call TriangleMeshParameters::disable_boundary_refinement()\n"
322  << "or enable internal boundary refinement (the default)\n";
323 
324  throw OomphLibError(error_stream.str().c_str(),
325  OOMPH_CURRENT_FUNCTION,
326  OOMPH_EXCEPTION_LOCATION);
327  }
328 
329  this->generic_constructor(
330  outer_boundary_polygon_pt,
331  internal_polygon_pt,
332  internal_open_curve_poly_pt,
333  element_area,
334  extra_holes_coordinates,
335  regions,
336  triangle_mesh_parameters.target_area_for_region(),
337  time_stepper_pt,
338  use_attributes,
339  refine_boundary,
340  refine_internal_boundary);
341 
342 #ifdef OOMPH_HAS_MPI
343 
344  // Before calling setup boundary coordinates check if the mesh is
345  // marked as distrbuted
346  if (triangle_mesh_parameters.is_mesh_distributed())
347  {
348  // Set the mesh as distributed by passing the communicator
349  this->set_communicator_pt(triangle_mesh_parameters.communicator_pt());
350  }
351 
352 #endif
353 
354  // Setup boundary coordinates for boundaries
355  unsigned nb = nboundary();
356  for (unsigned b = 0; b < nb; b++)
357  {
358  this->template setup_boundary_coordinates<ELEMENT>(b);
359  }
360 
361  // Snap it!
363  }
364 
365  /// A general-purpose construction function that builds the
366  /// mesh once the different specific constructors have assembled the
367  /// appropriate information.
369  Vector<TriangleMeshPolygon*>& outer_boundary_pt,
370  Vector<TriangleMeshPolygon*>& internal_polygon_pt,
371  Vector<TriangleMeshOpenCurve*>& open_polylines_pt,
372  const double& element_area,
373  Vector<Vector<double>>& extra_holes_coordinates,
374  std::map<unsigned, Vector<double>>& regions_coordinates,
375  std::map<unsigned, double>& regions_areas,
376  TimeStepper* time_stepper_pt,
377  const bool& use_attributes,
378  const bool& refine_boundary,
379  const bool& refine_internal_boundary)
380  {
381 #ifdef PARANOID
382 
383  if (element_area < 10e-14)
384  {
385  std::ostringstream warning_message;
386  warning_message
387  << "The current elements area was stated to (" << element_area
388  << ").\nThe current precision to generate the input to triangle "
389  << "is fixed to 14 digits\n\n";
390  OomphLibWarning(warning_message.str(),
391  OOMPH_CURRENT_FUNCTION,
392  OOMPH_EXCEPTION_LOCATION);
393  }
394 
395 #endif
396 
397  // Store the attribute flag
398  this->Use_attributes = use_attributes;
399 
400  // Store the timestepper used to build elements
401  this->Time_stepper_pt = time_stepper_pt;
402 
403  // Store outer polygon
404  this->Outer_boundary_pt = outer_boundary_pt;
405 
406  // Store internal polygons by copy constructor
407  this->Internal_polygon_pt = internal_polygon_pt;
408 
409  // Store internal polylines by copy constructor
410  this->Internal_open_curve_pt = open_polylines_pt;
411 
412  // Store the extra holes coordinates
413  this->Extra_holes_coordinates = extra_holes_coordinates;
414 
415  // Store the extra regions coordinates
416  this->Regions_coordinates = regions_coordinates;
417 
418  // Create the data structures required to call the triangulate function
419  TriangulateIO triangulate_io;
420  TriangulateIO triangulate_out;
421 
422  // Initialize TriangulateIO structure
424 
425  // Convert TriangleMeshPolyLine and TriangleMeshClosedCurvePolyLine
426  // to a triangulateio object
428  outer_boundary_pt,
429  internal_polygon_pt,
430  open_polylines_pt,
431  extra_holes_coordinates,
432  regions_coordinates,
433  regions_areas,
434  triangulate_io);
435 
436  // Initialize TriangulateIO structure
438 
439  // Input string for triangle
440  std::stringstream input_string_stream;
441  input_string_stream.precision(14);
442  input_string_stream.setf(std::ios_base::fixed, std::ios_base::floatfield);
443 
444  // MH: Used to be:
445  // input_string_stream<<"-pA -a" << element_area << " -q30" << std::fixed;
446  // The repeated -a allows the specification of areas for different
447  // regions (if any)
448  input_string_stream << "-pA -a -a" << element_area << " -q30"
449  << std::fixed;
450 
451  // Verify if creation of new points on boundaries is allowed
453  {
454  input_string_stream << " -YY";
455  }
456 
457  // Suppress insertion of additional points on outer boundary
458  if (refine_boundary == false)
459  {
460  input_string_stream << "-Y";
461 
462  // Add the extra flag to suppress additional points on interior segments
463  if (refine_internal_boundary == false)
464  {
465  input_string_stream << "Y";
466  }
467  }
468 
469  // Convert the Input string in *char required by the triangulate function
470  char triswitches[100];
471  sprintf(triswitches, "%s", input_string_stream.str().c_str());
472 
473  // Build the mesh using triangulate function
474  triangulate(triswitches, &triangulate_io, &triangulate_out, 0);
475 
476 #ifdef OOMPH_HAS_FPUCONTROLH
477  // Reset flags that are tweaked by triangle; can cause nasty crashes
478  fpu_control_t cw = (_FPU_DEFAULT & ~_FPU_EXTENDED) | _FPU_DOUBLE;
479  _FPU_SETCW(cw);
480 #endif
481 
482  // Build scaffold
483  TriangleScaffoldMesh* tmp_mesh_pt =
484  new TriangleScaffoldMesh(triangulate_out);
485 
486  // If we have filled holes then we must use the attributes
487  if (!regions_coordinates.empty())
488  {
489  // Convert mesh from scaffold to actual mesh
490  build_from_scaffold(tmp_mesh_pt, time_stepper_pt, true);
491 
492  // Record the attribute flag
493  this->Use_attributes = true;
494  }
495  // Otherwise use what was asked
496  else
497  {
498  // Convert mesh from scaffold to actual mesh
499  build_from_scaffold(tmp_mesh_pt, time_stepper_pt, use_attributes);
500  }
501 
502  // Kill the scaffold
503  delete tmp_mesh_pt;
504  tmp_mesh_pt = 0;
505 
506  // Cleanup but leave hole and regions alone since it's still used
507  bool clear_hole_data = false;
508  TriangleHelper::clear_triangulateio(triangulate_io, clear_hole_data);
509  TriangleHelper::clear_triangulateio(triangulate_out, clear_hole_data);
510  }
511 
512 #endif // OOMPH_HAS_TRIANGLE_LIB
513 
514  /// Broken copy constructor
516 
517  /// Broken assignment operator
518  void operator=(const QuadFromTriangleMesh&) = delete;
519 
520  /// Empty destructor
522  {
523 #ifdef OOMPH_HAS_TRIANGLE_LIB
524 
525  std::set<TriangleMeshCurveSection*>::iterator it_polyline;
526  for (it_polyline = Free_curve_section_pt.begin();
527  it_polyline != Free_curve_section_pt.end();
528  it_polyline++)
529  {
530  delete (*it_polyline);
531  }
532 
533  std::set<TriangleMeshPolygon*>::iterator it_polygon;
534  for (it_polygon = Free_polygon_pt.begin();
535  it_polygon != Free_polygon_pt.end();
536  it_polygon++)
537  {
538  delete (*it_polygon);
539  }
540 
541  std::set<TriangleMeshOpenCurve*>::iterator it_open_polyline;
542  for (it_open_polyline = Free_open_curve_pt.begin();
543  it_open_polyline != Free_open_curve_pt.end();
544  it_open_polyline++)
545  {
546  delete (*it_open_polyline);
547  }
548 
549 #endif
550  }
551 
552  /// Build the quad mesh from the given scaffold mesh
553  void build_from_scaffold(TriangleScaffoldMesh* tmp_mesh_pt,
554  TimeStepper* time_stepper_pt,
555  const bool& use_attributes);
556 
557  /// Timestepper used to build elements
559 
560  /// Boolean flag to indicate whether to use attributes or not (required
561  /// for multidomain meshes)
563  };
564 
565 
566  /// /////////////////////////////////////////////////////////////////
567  /// /////////////////////////////////////////////////////////////////
568  /// /////////////////////////////////////////////////////////////////
569 
570 
571  //=========================================================================
572  /// Unstructured refineable QuadFromTriangleMesh
573  //=========================================================================
574  template<class ELEMENT>
576  : public virtual QuadFromTriangleMesh<ELEMENT>,
577  public virtual RefineableQuadMesh<ELEMENT>
578  {
579  public:
580 #ifdef OOMPH_HAS_TRIANGLE_LIB
581 
582  /// Build mesh, based on the specifications on
583  /// TriangleMeshParameters
585  TriangleMeshParameters& triangle_mesh_parameters,
586  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper)
587  : QuadFromTriangleMesh<ELEMENT>(triangle_mesh_parameters, time_stepper_pt)
588  {
589  this->setup_quadtree_forest();
590  }
591 
592 #endif
593 
594  /// Refine mesh uniformly
595  virtual void refine_uniformly()
596  {
597  DocInfo doc_info;
598  doc_info.directory() = "";
599  doc_info.disable_doc();
600  refine_uniformly(doc_info);
601  }
602 
603  /// Refine mesh uniformly and doc process
604  void refine_uniformly(DocInfo& doc_info)
605  {
606  // Find the number of elements in the mesh
607  unsigned nelem = this->nelement();
608 
609  // Set the element error to something big
610  Vector<double> elem_error(nelem, DBL_MAX);
611 
612  // Refine everything
613  adapt(elem_error);
614  }
615 
616  /// Overload the adapt function (to ensure nodes are snapped to the
617  /// boundary)
618  void adapt(const Vector<double>& elem_error);
619 
620  /// Build mesh, based on the polyfiles
622  const std::string& node_file_name,
623  const std::string& element_file_name,
624  const std::string& poly_file_name,
625  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper)
626  : QuadFromTriangleMesh<ELEMENT>(
627  node_file_name, element_file_name, poly_file_name, time_stepper_pt)
628  {
629  this->setup_quadtree_forest();
630  }
631 
632  /// Empty Destructor
634  };
635 
636 
637  /// /////////////////////////////////////////////////////////////////
638  /// /////////////////////////////////////////////////////////////////
639  /// /////////////////////////////////////////////////////////////////
640 
641 
642  //=========================================================================
643  /// Unstructured QuadFromTriangleMesh upgraded to solid mesh
644  //=========================================================================
645  template<class ELEMENT>
647  : public virtual QuadFromTriangleMesh<ELEMENT>,
648  public virtual SolidMesh
649  {
650  public:
652  const std::string& node_file_name,
653  const std::string& element_file_name,
654  const std::string& poly_file_name,
655  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper,
656  const bool& use_attributes = false)
657  : QuadFromTriangleMesh<ELEMENT>(node_file_name,
658  element_file_name,
659  poly_file_name,
660  time_stepper_pt,
661  use_attributes)
662  {
663  // Assign the Lagrangian coordinates
664  set_lagrangian_nodal_coordinates();
665  }
666 
667 #ifdef OOMPH_HAS_TRIANGLE_LIB
668 
669  /// Build mesh, based on closed curve that specifies
670  /// the outer boundary of the domain and any number of internal
671  /// clsed curves. Specify target area for uniform element size.
673  TriangleMeshParameters& triangle_mesh_parameters,
674  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper)
675  : QuadFromTriangleMesh<ELEMENT>(triangle_mesh_parameters, time_stepper_pt)
676  {
677  // Assign the Lagrangian coordinates
678  set_lagrangian_nodal_coordinates();
679  }
680 
681 #endif
682 
683  /// Empty Destructor
685  };
686 
687 
688  /// /////////////////////////////////////////////////////////////////////
689  /// /////////////////////////////////////////////////////////////////////
690  /// /////////////////////////////////////////////////////////////////////
691 
692 
693  //=========================================================================
694  /// Unstructured refineable QuadFromTriangleMesh upgraded to solid mesh
695  //=========================================================================
696  template<class ELEMENT>
699  public virtual SolidMesh
700  {
701  public:
702  /// Build mesh from specified triangulation and associated
703  /// target areas for elements in it.
705  const std::string& node_file_name,
706  const std::string& element_file_name,
707  const std::string& poly_file_name,
708  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper,
709  const bool& use_attributes = false)
710  : RefineableQuadFromTriangleMesh<ELEMENT>(node_file_name,
711  element_file_name,
712  poly_file_name,
713  time_stepper_pt,
714  use_attributes)
715  {
716  // Assign the Lagrangian coordinates
717  set_lagrangian_nodal_coordinates();
718  }
719 
720 #ifdef OOMPH_HAS_TRIANGLE_LIB
721 
722  /// Build mesh, based on the specifications on
723  /// TriangleMeshParameter
725  TriangleMeshParameters& triangle_mesh_parameters,
726  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper)
727  : QuadFromTriangleMesh<ELEMENT>(triangle_mesh_parameters,
728  time_stepper_pt),
729  RefineableQuadFromTriangleMesh<ELEMENT>(triangle_mesh_parameters,
730  time_stepper_pt)
731  {
732  // Assign the Lagrangian coordinates
733  set_lagrangian_nodal_coordinates();
734  }
735 
736 #endif
737 
738  /// Empty Destructor
740  };
741 
742 
743 } // namespace oomph
744 
745 
746 #endif // OOMPH_QUAD_FROM_TRIANGLE_MESH_HEADER
e
Definition: cfortran.h:571
cstr elem_len * i
Definition: cfortran.h:603
Information for documentation of results: Directory and file number to enable output in the form RESL...
void disable_doc()
Disable documentation.
std::string directory() const
Output directory.
static Steady< 0 > Default_TimeStepper
Default Steady Timestepper, to be used in default arguments to Mesh constructors.
Definition: mesh.h:75
unsigned nboundary() const
Return number of boundaries.
Definition: mesh.h:827
void set_communicator_pt(OomphCommunicator *comm_pt)
Function to set communicator (mesh is assumed to be distributed if the communicator pointer is non-nu...
Definition: mesh.h:1615
An OomphLibError object which should be thrown when an run-time error is encountered....
An OomphLibWarning object which should be created as a temporary object to issue a warning....
Quad mesh built on top of triangle scaffold mesh coming from the triangle mesh generator Triangle....
void generic_constructor(Vector< TriangleMeshPolygon * > &outer_boundary_pt, Vector< TriangleMeshPolygon * > &internal_polygon_pt, Vector< TriangleMeshOpenCurve * > &open_polylines_pt, const double &element_area, Vector< Vector< double >> &extra_holes_coordinates, std::map< unsigned, Vector< double >> &regions_coordinates, std::map< unsigned, double > &regions_areas, TimeStepper *time_stepper_pt, const bool &use_attributes, const bool &refine_boundary, const bool &refine_internal_boundary)
A general-purpose construction function that builds the mesh once the different specific constructors...
TimeStepper * Time_stepper_pt
Timestepper used to build elements.
QuadFromTriangleMesh(TriangleMeshParameters &triangle_mesh_parameters, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Build mesh, based on the specifications on TriangleMeshParameters. All the actual work is done in Uns...
void operator=(const QuadFromTriangleMesh &)=delete
Broken assignment operator.
QuadFromTriangleMesh(const std::string &node_file_name, const std::string &element_file_name, const std::string &poly_file_name, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper, const bool &use_attributes=false, const bool &allow_automatic_creation_of_vertices_on_boundaries=true)
Constructor with the input files.
bool Use_attributes
Boolean flag to indicate whether to use attributes or not (required for multidomain meshes)
void build_from_scaffold(TriangleScaffoldMesh *tmp_mesh_pt, TimeStepper *time_stepper_pt, const bool &use_attributes)
Build the quad mesh from the given scaffold mesh.
QuadFromTriangleMesh(const QuadFromTriangleMesh &dummy)=delete
Broken copy constructor.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: quad_mesh.h:57
///////////////////////////////////////////////////////////////// ///////////////////////////////////...
virtual void refine_uniformly()
Refine mesh uniformly.
RefineableQuadFromTriangleMesh(TriangleMeshParameters &triangle_mesh_parameters, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Build mesh, based on the specifications on TriangleMeshParameters.
RefineableQuadFromTriangleMesh(const std::string &node_file_name, const std::string &element_file_name, const std::string &poly_file_name, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Build mesh, based on the polyfiles.
void refine_uniformly(DocInfo &doc_info)
Refine mesh uniformly and doc process.
Intermediate mesh class that implements the mesh adaptation functions specified in the TreeBasedRefin...
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
RefineableSolidQuadFromTriangleMesh(TriangleMeshParameters &triangle_mesh_parameters, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Build mesh, based on the specifications on TriangleMeshParameter.
RefineableSolidQuadFromTriangleMesh(const std::string &node_file_name, const std::string &element_file_name, const std::string &poly_file_name, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper, const bool &use_attributes=false)
Build mesh from specified triangulation and associated target areas for elements in it.
General SolidMesh class.
Definition: mesh.h:2562
///////////////////////////////////////////////////////////////// ///////////////////////////////////...
SolidQuadFromTriangleMesh(TriangleMeshParameters &triangle_mesh_parameters, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Build mesh, based on closed curve that specifies the outer boundary of the domain and any number of i...
SolidQuadFromTriangleMesh(const std::string &node_file_name, const std::string &element_file_name, const std::string &poly_file_name, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper, const bool &use_attributes=false)
////////////////////////////////////////////////////////////////////// //////////////////////////////...
Definition: timesteppers.h:231
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
Vector< TriangleMeshClosedCurve * > internal_closed_curve_pt() const
Helper function for getting the internal closed boundaries.
bool is_automatic_creation_of_vertices_on_boundaries_allowed()
Returns the status of the variable Allow_automatic_creation_of_vertices_on_boundaries.
Vector< Vector< double > > extra_holes_coordinates() const
Helper function for getting the extra holes.
std::map< unsigned, double > & target_area_for_region()
Helper function for getting access to the region's target areas.
bool is_mesh_distributed() const
Boolean to indicate if Mesh has been distributed.
bool is_use_attributes() const
Helper function for getting the status of use_attributes variable.
Vector< TriangleMeshClosedCurve * > outer_boundary_pt() const
Helper function for getting the outer boundary.
OomphCommunicator * communicator_pt() const
Read-only access fct to communicator (Null if mesh is not distributed)
double element_area() const
Helper function for getting the element area.
Vector< TriangleMeshOpenCurve * > internal_open_curves_pt() const
Helper function for getting the internal open boundaries.
std::map< unsigned, Vector< double > > & regions_coordinates()
Helper function for getting access to the regions coordinates.
bool is_boundary_refinement_allowed() const
Helper function for getting the status of boundary refinement.
bool is_internal_boundary_refinement_allowed() const
Helper function for getting the status of boundary refinement.
Triangle Mesh that is based on input files generated by the triangle mesh generator Triangle.
Contains functions which define the geometry of the mesh, i.e. regions, boundaries,...
Vector< TriangleMeshOpenCurve * > Internal_open_curve_pt
Vector of open polylines that define internal curves.
bool is_automatic_creation_of_vertices_on_boundaries_allowed()
Returns the status of the variable Allow_automatic_creation_of_vertices_on_boundaries.
Vector< TriangleMeshPolygon * > Outer_boundary_pt
Polygon that defines outer boundaries.
Vector< TriangleMeshPolygon * > Internal_polygon_pt
Vector of polygons that define internal polygons.
void snap_nodes_onto_geometric_objects()
Snap the boundary nodes onto any curvilinear boundaries defined by geometric objects.
TriangleMeshOpenCurve * create_open_curve_with_polyline_helper(TriangleMeshOpenCurve *open_curve_pt, unsigned &max_bnd_id_local)
Helper function that creates and returns an open curve with the polyline representation of its consti...
std::map< unsigned, Vector< double > > Regions_coordinates
Storage for extra coordinates for regions. The key on the map is the region id.
Vector< Vector< double > > Extra_holes_coordinates
Storage for extra coordinates for holes.
bool Allow_automatic_creation_of_vertices_on_boundaries
Flag to indicate whether the automatic creation of vertices along the boundaries by Triangle is allow...
std::set< TriangleMeshOpenCurve * > Free_open_curve_pt
A set that contains the open curves created by this object therefore it is necessary to free their as...
TriangleMeshPolygon * closed_curve_to_polygon_helper(TriangleMeshClosedCurve *closed_curve_pt, unsigned &max_bnd_id_local)
Helper function that returns a polygon representation for the given closed curve, it also computes th...
void build_triangulateio(Vector< TriangleMeshPolygon * > &outer_polygons_pt, Vector< TriangleMeshPolygon * > &internal_polygons_pt, Vector< TriangleMeshOpenCurve * > &open_curves_pt, Vector< Vector< double >> &extra_holes_coordinates, std::map< unsigned, Vector< double >> &regions_coordinates, std::map< unsigned, double > &regions_areas, TriangulateIO &triangulate_io)
Create TriangulateIO object from outer boundaries, internal boundaries, and open curves....
std::set< TriangleMeshCurveSection * > Free_curve_section_pt
A set that contains the curve sections created by this object therefore it is necessary to free their...
std::set< TriangleMeshPolygon * > Free_polygon_pt
A set that contains the polygons created by this object therefore it is necessary to free their assoc...
void set_geom_objects_and_coordinate_limits_for_open_curve(TriangleMeshOpenCurve *input_open_curve_pt)
Stores the geometric objects associated to the curve sections that compound the open curve....
void set_geom_objects_and_coordinate_limits_for_close_curve(TriangleMeshClosedCurve *input_closed_curve_pt)
Stores the geometric objects associated to the curve sections that compound the closed curve....
A slight extension to the standard template vector class so that we can include "graceful" array rang...
Definition: Vector.h:58
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
void initialise_triangulateio(TriangulateIO &triangle_io)
Initialise TriangulateIO structure.
void clear_triangulateio(TriangulateIO &triangulate_io, const bool &clear_hole_data)
Clear TriangulateIO structure.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void triangulate(char *triswitches, struct oomph::TriangulateIO *in, struct oomph::TriangulateIO *out, struct oomph::TriangulateIO *vorout)
The Triangle data structure, modified from the triangle.h header supplied with triangle 1....