OOFEM 3.0
Loading...
Searching...
No Matches
edge2d.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#ifndef edge2d_h
36#define edge2d_h
37
38#include "intarray.h"
39
40namespace oofem {
42class IntArray;
43
49class Edge2D
50{
51private:
53 std :: pair< int, int >nodeNumbers;
54
55public:
57 Edge2D(int node1, int node2);
59 virtual ~Edge2D();
61 int giveFirstNodeNumber() { return nodeNumbers.first; }
63 int giveSecondNodeNumber() { return nodeNumbers.second; }
65 virtual bool operator==(const Edge2D &right);
66};
67
77class AlphaEdge2D : public Edge2D
78{
79public:
81 AlphaEdge2D(int node1, int node2, double _length);
83 virtual ~AlphaEdge2D();
84
86 void setOuterAlphaBound(double alphaMin) { outerAlphaBound = alphaMin; }
88 void setInnerAlphaBound(double alphaMax) { innerAlphaBound = alphaMax; }
94 void setHullFlag(bool flag) { isOnConvexHull = flag; }
96 bool giveHullFlag() { return isOnConvexHull; }
98 void setSharing(int n, DelaunayTriangle *pTE);
100 DelaunayTriangle *giveShared(int n) { return sharedByTriangles [ n - 1 ]; }
102 double giveLength() { return length; }
103
104private:
114 double length;
115};
116}
117#endif
double giveOuterAlphaBound()
Returns the outer limit.
Definition edge2d.h:90
void setOuterAlphaBound(double alphaMin)
Sets the outer limit.
Definition edge2d.h:86
void setSharing(int n, DelaunayTriangle *pTE)
Stores DelaunayTriangle sharing the receiver.
Definition edge2d.C:70
DelaunayTriangle * giveShared(int n)
Returns DelaunayTriangle sharing receiver.
Definition edge2d.h:100
void setHullFlag(bool flag)
Sets the convex hull property.
Definition edge2d.h:94
AlphaEdge2D(int node1, int node2, double _length)
Constructor.
Definition edge2d.C:55
bool giveHullFlag()
Returns true if the edge lies on convex hull, false otherwise.
Definition edge2d.h:96
DelaunayTriangle * sharedByTriangles[2]
Triangles which share the alphaEdge.
Definition edge2d.h:112
double outerAlphaBound
Bottom (outer) limit for alpha shape, for smaller values of alpha the edge lies outside.
Definition edge2d.h:108
double length
Length of edge is stored in order to allow variable alpha value.
Definition edge2d.h:114
bool isOnConvexHull
Convex hull flag means edge is not shared by two triangle.
Definition edge2d.h:106
double innerAlphaBound
Top (inner) limit for alpha shape, for greater values of alpha the edge lies inside of the shape.
Definition edge2d.h:110
double giveLength()
Returns length of the receiver.
Definition edge2d.h:102
void setInnerAlphaBound(double alphaMax)
Sets the inner limit.
Definition edge2d.h:88
virtual ~AlphaEdge2D()
Destructor.
Definition edge2d.C:66
double giveInnerAlphaBound()
Returns the inner limit.
Definition edge2d.h:92
Edge2D(int node1, int node2)
Constructor.
Definition edge2d.C:39
virtual ~Edge2D()
Destructor.
Definition edge2d.C:43
std ::pair< int, int > nodeNumbers
Global node numbers.
Definition edge2d.h:53
int giveSecondNodeNumber()
Gives the number of the second node.
Definition edge2d.h:63
virtual bool operator==(const Edge2D &right)
Compares receiver with passed Edge2D. Returns true if node numbers are equal, false otherwise.
Definition edge2d.C:47
int giveFirstNodeNumber()
Gives the number of the first node.
Definition edge2d.h:61

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