OOFEM 3.0
Loading...
Searching...
No Matches
delaunaytriangulator.h
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34
35// *******************************
36// *** DELAUNAY MESH GENERATOR ***
37// *******************************
38
39// based on core/mesherinterface.h
40
41
42#ifndef delaunaytrinagulator_h
43#define delaunaytrinagulator_h
44
45#include <list>
46#include "contextioresulttype.h"
47#include "timer.h"
48#include "octreelocalizert.h"
50#include "intarray.h"
51#include "edge2d.h"
52
53namespace oofem {
54class Domain;
55class TimeStep;
56class IntArray;
57class FloatArray;
58class InsertionData;
59
67{
68protected:
71 // Value of alpha for the boundary recognition via alpha shape algorithm
72 double alphaValue;
73 // Number of nodes
74 int nnode;
75 // Option 2: setting bounds for computed Alpha
76 //double minAlpha;
77 //double maxAlpha;
78
88 std :: list< DelaunayTriangle * >generalTriangleList;
89
91 std :: list< AlphaEdge2D * >alphaShapeEdgeList;
92
94 std :: list< AlphaEdge2D * >edgeList;
95
98
99public:
101 DelaunayTriangulator(Domain *d, double setAlpha);
104
106 void generateMesh();
107
108private:
110 void addUniqueEdgeToPolygon(Edge2D edge, std :: list< Edge2D > &polygon);
111
115 void writeMesh();
116
118 void computeAlphaComplex();
119
124
126 void giveAlphaShape();
127
129 void initializeTimers();
130
132 void findNonDelaunayTriangles(int insertedNode, InsertTriangleBasedOnCircumcircle &tInsert, std :: list< Edge2D > &polygon);
134 void meshPolygon(int insertedNode, InsertTriangleBasedOnCircumcircle &tInsert, std :: list< Edge2D > &polygon);
136 void giveTimeReport();
138 void cleanUpTriangleList();
139
142
143};
144} // end namespace oofem
145#endif // delaunaytrinagulator_h
Timer searchingTimer
Measures time needed by searching for non-delaunay triangles.
void addUniqueEdgeToPolygon(Edge2D edge, std ::list< Edge2D > &polygon)
Edge is added to the polygon only if it's not contained. Otherwise both are removed (edge shared by t...
void giveAlphaShape()
Iterates through the edgeList container and compares alpha-value with alphaEdge bounds....
Timer meshingTimer
Measures overall time of triangulation procedure.
void meshPolygon(int insertedNode, InsertTriangleBasedOnCircumcircle &tInsert, std ::list< Edge2D > &polygon)
Retriangulates the polygon.
std ::list< AlphaEdge2D * > alphaShapeEdgeList
Contains resulting alpha-shape in form of a list.
void computeAlphaComplex()
Reads the triangulation and fills tha edgeList container with alpha-shape edges and set their bounds.
void giveTimeReport()
Prints the time report.
Timer polygonTimer
Measures time needed for identifying polygon to be retriangulated.
AlphaEdge2D * giveBackEdgeIfAlreadyContainedInList(AlphaEdge2D &alphaEdge)
void computeBBXBasedOnNodeData(BoundingBox &BBX)
Calculates the bounding box base on the domain's nodes.
std ::list< AlphaEdge2D * > edgeList
contains all edges of the triangulation
OctreeSpatialLocalizerT< DelaunayTriangle * > triangleOctree
Octree with Delaunay triangles allowing fast search.
Domain * domain
Domain of the PFEM problem containing nodes to be triangulated.
void buildInitialBBXMesh(InsertTriangleBasedOnCircumcircle &tInsert)
Identifies the bounding box of pfemparticles and creates initial triangulation consisting of 2 triang...
DelaunayTriangulator(Domain *d, double setAlpha)
Constructor.
void cleanUpTriangleList()
Iterates through generalTringleList und removes non-valid ones or those containing bounding box nodes...
void writeMesh()
Writes the mesh into the domain by creating new tr1_2d_pfem elements and prescribes zero-pressure bou...
void initializeTimers()
Initializes Timers and.
Timer creativeTimer
Measures time needed for creating new Delaunay triangles.
void findNonDelaunayTriangles(int insertedNode, InsertTriangleBasedOnCircumcircle &tInsert, std ::list< Edge2D > &polygon)
Looks for non-Delaunay triangles in octree and creates a polygon.
std ::list< DelaunayTriangle * > generalTriangleList
Contains all triangles (even not valid).

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011