quad_from_triangle_mesh.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-2025 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_HEADER
28#define OOMPH_QUAD_FROM_TRIANGLE_MESH_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.h"
54#include "generic/Qelements.h"
55
56
57////////////////////////////////////////////////////////////////////
58////////////////////////////////////////////////////////////////////
59////////////////////////////////////////////////////////////////////
60
61
62namespace 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,
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
102
103 // Store Timestepper used to build elements
104 this->Time_stepper_pt = time_stepper_pt;
105
106 // Store the attributes
108
109 // Build scaffold
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 {
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
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(),
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
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
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)
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
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
233 triangle_mesh_parameters.internal_open_curves_pt();
234
235 // Find the number of internal open curves
237
238 // Create the storage for the polylines that define the open boundaries
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)
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;
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(),
327 }
328
330 outer_boundary_polygon_pt,
333 element_area,
334 extra_holes_coordinates,
335 regions,
336 triangle_mesh_parameters.target_area_for_region(),
337 time_stepper_pt,
339 refine_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,
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;
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";
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
421
422 // Initialize TriangulateIO structure
424
425 // Convert TriangleMeshPolyLine and TriangleMeshClosedCurvePolyLine
426 // to a triangulateio object
428 outer_boundary_pt,
431 extra_holes_coordinates,
432 regions_coordinates,
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];
472 sizeof(triswitches),
473 "%s",
474 input_string_stream.str().c_str());
475
476 // Build the mesh using triangulate function
478
479#ifdef OOMPH_HAS_FPUCONTROLH
480 // Reset flags that are tweaked by triangle; can cause nasty crashes
482 _FPU_SETCW(cw);
483#endif
484
485 // Build scaffold
488
489 // If we have filled holes then we must use the attributes
490 if (!regions_coordinates.empty())
491 {
492 // Convert mesh from scaffold to actual mesh
493 build_from_scaffold(tmp_mesh_pt, time_stepper_pt, true);
494
495 // Record the attribute flag
496 this->Use_attributes = true;
497 }
498 // Otherwise use what was asked
499 else
500 {
501 // Convert mesh from scaffold to actual mesh
502 build_from_scaffold(tmp_mesh_pt, time_stepper_pt, use_attributes);
503 }
504
505 // Kill the scaffold
506 delete tmp_mesh_pt;
507 tmp_mesh_pt = 0;
508
509 // Cleanup but leave hole and regions alone since it's still used
510 bool clear_hole_data = false;
513 }
514
515#endif // OOMPH_HAS_TRIANGLE_LIB
516
517 /// Broken copy constructor
519
520 /// Broken assignment operator
521 void operator=(const QuadFromTriangleMesh&) = delete;
522
523 /// Empty destructor
525 {
526#ifdef OOMPH_HAS_TRIANGLE_LIB
527
528 std::set<TriangleMeshCurveSection*>::iterator it_polyline;
529 for (it_polyline = Free_curve_section_pt.begin();
531 it_polyline++)
532 {
533 delete (*it_polyline);
534 }
535
536 std::set<TriangleMeshPolygon*>::iterator it_polygon;
537 for (it_polygon = Free_polygon_pt.begin();
539 it_polygon++)
540 {
541 delete (*it_polygon);
542 }
543
544 std::set<TriangleMeshOpenCurve*>::iterator it_open_polyline;
548 {
549 delete (*it_open_polyline);
550 }
551
552#endif
553 }
554
555 /// Build the quad mesh from the given scaffold mesh
557 TimeStepper* time_stepper_pt,
558 const bool& use_attributes);
559
560 /// Timestepper used to build elements
562
563 /// Boolean flag to indicate whether to use attributes or not (required
564 /// for multidomain meshes)
566 };
567
568
569 ////////////////////////////////////////////////////////////////////
570 ////////////////////////////////////////////////////////////////////
571 ////////////////////////////////////////////////////////////////////
572
573
574 //=========================================================================
575 /// Unstructured refineable QuadFromTriangleMesh
576 //=========================================================================
577 template<class ELEMENT>
579 : public virtual QuadFromTriangleMesh<ELEMENT>,
580 public virtual RefineableQuadMesh<ELEMENT>
581 {
582 public:
583#ifdef OOMPH_HAS_TRIANGLE_LIB
584
585 /// Build mesh, based on the specifications on
586 /// TriangleMeshParameters
589 TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper)
590 : QuadFromTriangleMesh<ELEMENT>(triangle_mesh_parameters, time_stepper_pt)
591 {
592 this->setup_quadtree_forest();
593 }
594
595#endif
596
597 /// Refine mesh uniformly
598 virtual void refine_uniformly()
599 {
600 DocInfo doc_info;
601 doc_info.directory() = "";
602 doc_info.disable_doc();
603 refine_uniformly(doc_info);
604 }
605
606 /// Refine mesh uniformly and doc process
607 void refine_uniformly(DocInfo& doc_info)
608 {
609 // Find the number of elements in the mesh
610 unsigned nelem = this->nelement();
611
612 // Set the element error to something big
614
615 // Refine everything
616 adapt(elem_error);
617 }
618
619 /// Overload the adapt function (to ensure nodes are snapped to the
620 /// boundary)
621 void adapt(const Vector<double>& elem_error);
622
623 /// Build mesh, based on the polyfiles
625 const std::string& node_file_name,
626 const std::string& element_file_name,
627 const std::string& poly_file_name,
628 TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper)
629 : QuadFromTriangleMesh<ELEMENT>(
631 {
632 this->setup_quadtree_forest();
633 }
634
635 /// Empty Destructor
637 };
638
639
640 ////////////////////////////////////////////////////////////////////
641 ////////////////////////////////////////////////////////////////////
642 ////////////////////////////////////////////////////////////////////
643
644
645 //=========================================================================
646 /// Unstructured QuadFromTriangleMesh upgraded to solid mesh
647 //=========================================================================
648 template<class ELEMENT>
650 : public virtual QuadFromTriangleMesh<ELEMENT>,
651 public virtual SolidMesh
652 {
653 public:
655 const std::string& node_file_name,
656 const std::string& element_file_name,
657 const std::string& poly_file_name,
658 TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper,
659 const bool& use_attributes = false)
663 time_stepper_pt,
665 {
666 // Assign the Lagrangian coordinates
667 set_lagrangian_nodal_coordinates();
668 }
669
670#ifdef OOMPH_HAS_TRIANGLE_LIB
671
672 /// Build mesh, based on closed curve that specifies
673 /// the outer boundary of the domain and any number of internal
674 /// clsed curves. Specify target area for uniform element size.
677 TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper)
678 : QuadFromTriangleMesh<ELEMENT>(triangle_mesh_parameters, time_stepper_pt)
679 {
680 // Assign the Lagrangian coordinates
681 set_lagrangian_nodal_coordinates();
682 }
683
684#endif
685
686 /// Empty Destructor
688 };
689
690
691 ////////////////////////////////////////////////////////////////////////
692 ////////////////////////////////////////////////////////////////////////
693 ////////////////////////////////////////////////////////////////////////
694
695
696 //=========================================================================
697 /// Unstructured refineable QuadFromTriangleMesh upgraded to solid mesh
698 //=========================================================================
699 template<class ELEMENT>
702 public virtual SolidMesh
703 {
704 public:
705 /// Build mesh from specified triangulation and associated
706 /// target areas for elements in it.
708 const std::string& node_file_name,
709 const std::string& element_file_name,
710 const std::string& poly_file_name,
711 TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper,
712 const bool& use_attributes = false)
716 time_stepper_pt,
718 {
719 // Assign the Lagrangian coordinates
720 set_lagrangian_nodal_coordinates();
721 }
722
723#ifdef OOMPH_HAS_TRIANGLE_LIB
724
725 /// Build mesh, based on the specifications on
726 /// TriangleMeshParameter
729 TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper)
731 time_stepper_pt),
733 time_stepper_pt)
734 {
735 // Assign the Lagrangian coordinates
736 set_lagrangian_nodal_coordinates();
737 }
738
739#endif
740
741 /// Empty Destructor
743 };
744
745
746} // namespace oomph
747
748
749#include "quad_from_triangle_mesh.template.cc" // OOMPH_QUAD_FROM_TRIANGLE_MESH_HEADER
750#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.
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
Definition elements.cc:4320
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:835
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:1623
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....
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 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...
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.
Base class for quad meshes (meshes made of 2D quad elements).
Definition quad_mesh.h:57
Unstructured refineable QuadFromTriangleMesh.
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.
virtual ~RefineableQuadFromTriangleMesh()
Empty Destructor.
void refine_uniformly(DocInfo &doc_info)
Refine mesh uniformly and doc process.
Unstructured refineable QuadFromTriangleMesh upgraded to solid mesh.
virtual ~RefineableSolidQuadFromTriangleMesh()
Empty Destructor.
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:2570
Unstructured QuadFromTriangleMesh upgraded to solid mesh.
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)
virtual ~SolidQuadFromTriangleMesh()
Empty Destructor.
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
Base class for time-stepping schemes. Timestepper provides an approximation of the temporal derivativ...
Helper object for dealing with the parameters used for the TriangleMesh objects.
bool is_automatic_creation_of_vertices_on_boundaries_allowed()
Returns the status of the variable Allow_automatic_creation_of_vertices_on_boundaries.
OomphCommunicator * communicator_pt() const
Read-only access fct to communicator (Null if mesh is not distributed)
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...
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...
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...
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 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....
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
void initialise_triangulateio(TriangulateIO &triangle_io)
Initialise TriangulateIO structure.
void clear_triangulateio(TriangulateIO &triangulate_io, const bool &clear_hole_data)
Clear TriangulateIO structure.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
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....