faceAreaIntersectI.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2018 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30
32(
33 const face& f,
35)
36{
37 if (f.size() > 2)
38 {
39 const label v0 = 0;
40
41 face indices(3);
42
43 for (label i = 1; i < f.size() - 1; ++i)
44 {
45 indices[0] = f[v0];
46 indices[1] = f[i];
47 indices[2] = f[i + 1];
48 faces.append(indices);
49 }
50 }
51}
52
53
54inline void Foam::faceAreaIntersect::setTriPoints
55(
56 const point& a,
57 const point& b,
58 const point& c,
59 label& count,
61) const
62{
63 triPoints& tp = tris[count++];
64 tp[0] = a;
65 tp[1] = b;
66 tp[2] = c;
67}
68
69
70inline Foam::triPoints Foam::faceAreaIntersect::getTriPoints
71(
72 const pointField& points,
73 const face& f,
74 const bool reverse
75) const
76{
77 triPoints result;
78
79 if (reverse)
80 {
81 result[2] = points[f[0]];
82 result[1] = points[f[1]];
83 result[0] = points[f[2]];
84 }
85 else
86 {
87 result[0] = points[f[0]];
88 result[1] = points[f[1]];
89 result[2] = points[f[2]];
90 }
91
92 return result;
93}
94
95
96inline Foam::point Foam::faceAreaIntersect::planeIntersection
97(
98 const FixedList<scalar, 3>& d,
99 const triPoints& t,
100 const label negI,
101 const label posI
102) const
103{
104 scalar dp = d[posI];
105 scalar dn = d[negI];
106 return (dp*t[negI] - dn*t[posI])/(-dn + dp);
107}
108
109
110inline Foam::scalar Foam::faceAreaIntersect::triArea(const triPoints& t) const
111{
112 return mag(0.5*((t[1] - t[0])^(t[2] - t[0])));
113}
114
115
116inline Foam::vector Foam::faceAreaIntersect::triCentroid
117(
118 const triPoints& t
119) const
120{
121 return (t[0] + t[1] + t[2])/3;
122}
123
124
125// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
126
128{
129 return tol;
130}
131
132
134{
135 return cacheTriangulation_;
136}
137
138
141{
142 return triangles_;
143}
144
145
147{
148 return triangles_;
149}
150
151
152// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: FixedList.H:81
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static void triangleFan(const face &f, DynamicList< face > &faces)
Decompose face into triangle fan.
const DynamicList< triPoints > triangles() const
Const access to the triangulation.
bool cacheTriangulation() const
Const access to the cacheTriangulation flag.
static scalar & tolerance()
Fraction of local length scale to use as intersection tolerance.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Triangle storage. Null constructable (unfortunately triangle<point, point> is not)
Definition: triPoints.H:55
const pointField & points
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition: UListI.H:449
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList f(nPoints)
volScalarField & b
Definition: createFields.H:27