Loading...
Searching...
No Matches
tetIndicesI.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) 2021-2022 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// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30
33 celli_(-1),
34 facei_(-1),
35 tetPti_(-1)
36{}
37
38
39inline constexpr Foam::tetIndices::tetIndices
40(
41 label celli,
42 label facei,
43 label tetPointi
45:
46 celli_(celli),
47 facei_(facei),
48 tetPti_(tetPointi)
49{}
50
51
52// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
53
55(
56 const polyMesh& mesh,
57 const bool warn
58) const
59{
60 const Foam::face& f = mesh.faces()[facei_];
61
62 label faceBasePtI = mesh.tetBasePtIs()[facei_];
63
64 if (faceBasePtI < 0)
65 {
66 faceBasePtI = 0;
67
68 if (warn && nWarnings_ < maxNWarnings)
69 {
71 << "No base point for face " << facei_ << ", " << f
72 << ", produces a valid tet decomposition." << endl;
73
74 if (++nWarnings_ == maxNWarnings)
75 {
77 << "Suppressing further warnings." << endl;
78 }
79 }
80 }
81
82 const label facePtI = (tetPti_ + faceBasePtI) % f.size();
83 const label faceOtherPtI = f.fcIndex(facePtI);
84
85 if (mesh.faceOwner()[facei_] != celli_)
86 {
87 // Neighbour: return flipped face
88 return triFace(f[faceBasePtI], f[faceOtherPtI], f[facePtI]);
89 }
90
91 return triFace(f[faceBasePtI], f[facePtI], f[faceOtherPtI]);
92}
93
94
96(
97 const polyMesh& mesh,
98 const bool warn
99) const
100{
101 const Foam::face& f = mesh.faces()[facei_];
102
103 label faceBasePtI = mesh.tetBasePtIs()[facei_];
104
105 if (faceBasePtI < 0)
106 {
107 faceBasePtI = 0;
108
109 if (warn && nWarnings_ < maxNWarnings)
110 {
112 << "No base point for face " << facei_ << ", " << f
113 << ", produces a valid tet decomposition." << endl;
114
115 if (++nWarnings_ == maxNWarnings)
116 {
117 Warning
118 << "Suppressing further warnings." << endl;
119 }
120 }
121 }
122
123 const label facePtI = (tetPti_ + faceBasePtI) % f.size();
124 const label faceOtherPtI = f.fcIndex(facePtI);
125
126 if (mesh.faceOwner()[facei_] != celli_)
127 {
128 // Neighbour: return flipped face
129 return triFace(faceBasePtI, faceOtherPtI, facePtI);
130 }
131
132 return triFace(faceBasePtI, facePtI, faceOtherPtI);
133}
134
135
136inline Foam::tetPointRef Foam::tetIndices::tet(const polyMesh& mesh) const
137{
138 const pointField& pts = mesh.points();
139 const triFace tri = faceTriIs(mesh);
140
141 return tetPointRef
142 (
143 mesh.cellCentres()[celli_],
144 pts[tri[0]],
145 pts[tri[1]],
146 pts[tri[2]]
147 );
148}
149
150
152{
153 const pointField& pts = mesh.oldPoints();
154 const triFace tri = faceTriIs(mesh);
155
156 return tetPointRef
157 (
158 mesh.oldCellCentres()[celli_],
159 pts[tri[0]],
160 pts[tri[1]],
161 pts[tri[2]]
162 );
163}
164
165
167(
168 const polyMesh& mesh,
169 const barycentric& bary
170) const
171{
172 const pointField& pts = mesh.points();
173 const triFace tri = faceTriIs(mesh);
174
175 const point& a = mesh.cellCentres()[celli_];
176 const point& b = pts[tri[0]];
177 const point& c = pts[tri[1]];
178 const point& d = pts[tri[2]];
179
180 return point
181 (
182 (bary.a()*a.x() + bary.b()*b.x() + bary.c()*c.x() + bary.d()*d.x()),
183 (bary.a()*a.y() + bary.b()*b.y() + bary.c()*c.y() + bary.d()*d.y()),
184 (bary.a()*a.z() + bary.b()*b.z() + bary.c()*c.z() + bary.d()*d.z())
185 );
186}
187
190{
191 return faceTriIs(mesh).tri(mesh.points());
192}
193
194
196(
198) const
199{
200 return faceTriIs(mesh).tri(mesh.oldPoints());
201}
202
203
204inline Foam::label Foam::tetIndices::compare
205(
206 const tetIndices& a,
207 const tetIndices& b
208) noexcept
209{
210 label diff;
211 return
212 (
213 ((diff = (a.cell() - b.cell())) != 0) ? diff
214 : ((diff = (a.face() - b.face())) != 0) ? diff
215 : (a.tetPt() - b.tetPt())
216 );
217}
218
219
220// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
221
222inline bool Foam::operator==(const tetIndices& a, const tetIndices& b) noexcept
223{
224 // Possibly slightly faster version than compare
225 return
226 (
227 a.cell() == b.cell()
228 && a.face() == b.face()
229 && a.tetPt() == b.tetPt()
230 );
231}
232
233
234inline bool Foam::operator!=(const tetIndices& a, const tetIndices& b) noexcept
235{
236 return !(a == b);
237}
238
239
240// ************************************************************************* //
const Cmpt & b() const noexcept
const Cmpt & d() const noexcept
const Cmpt & a() const noexcept
const Cmpt & c() const noexcept
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition UListI.H:99
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition tetIndices.H:79
tetPointRef oldTet(const polyMesh &mesh) const
The tet geometry for this tet (using old positions), where point0 is the cell centre.
triFace triIs(const polyMesh &mesh, const bool warn=true) const
Local indices corresponding to the tri on the face for this tet. The normal of the tri points out of ...
Definition tetIndicesI.H:89
static label compare(const tetIndices &a, const tetIndices &b) noexcept
Compare tetIndices for equality. Compares cell, face, tetPt elements in order, stopping at the first ...
triFace faceTriIs(const polyMesh &mesh, const bool warn=true) const
The indices corresponding to the tri on the face for this tet. The normal of the tri points out of th...
Definition tetIndicesI.H:48
tetPointRef tet(const polyMesh &mesh) const
The tet geometry for this tet, where point0 is the cell centre.
triPointRef faceTri(const polyMesh &mesh) const
The triangle geometry for the face for this tet. The normal of the tri points out of the cell.
constexpr tetIndices() noexcept
Default construct, with invalid labels (-1).
Definition tetIndicesI.H:24
triPointRef oldFaceTri(const polyMesh &mesh) const
The triangle geometry for the face for this tet (using old positions).
point barycentricToPoint(const polyMesh &mesh, const barycentric &bary) const
The x/y/z position for given barycentric coordinates (where point0 is the cell centre).
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition triFace.H:68
dynamicFvMesh & mesh
#define WarningInFunction
Report a warning using Foam::Warning.
bool operator!=(const eddy &a, const eddy &b)
Definition eddy.H:297
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition tetrahedron.H:72
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition triad.C:373
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
vectorField pointField
pointField is a vectorField.
messageStream Warning
Warning stream (stdout output on master, null elsewhere), with additional 'FOAM Warning' header text.
face triFace(3)
labelList f(nPoints)
const pointField & pts
volScalarField & b