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
31inline constexpr Foam::tetIndices::tetIndices() noexcept
32:
33 celli_(-1),
34 facei_(-1),
35 tetPti_(-1)
36{}
37
38
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()[face()];
61
62 label faceBasePtI = mesh.tetBasePtIs()[face()];
63
64 if (faceBasePtI < 0)
65 {
66 faceBasePtI = 0;
67
68 if (warn && nWarnings_ < maxNWarnings)
69 {
70 ++nWarnings_;
72 << "No base point for face " << face() << ", " << f
73 << ", produces a valid tet decomposition." << endl;
74 if (nWarnings_ == maxNWarnings)
75 {
77 << "Suppressing any further warnings." << endl;
78 }
79 }
80 }
81
82 label facePtI = (tetPti_ + faceBasePtI) % f.size();
83 label faceOtherPtI = f.fcIndex(facePtI);
84
85 if (mesh.faceOwner()[face()] != cell())
86 {
87 std::swap(facePtI, faceOtherPtI);
88 }
89
90 return triFace(f[faceBasePtI], f[facePtI], f[faceOtherPtI]);
91}
92
93
95(
96 const polyMesh& mesh,
97 const bool warn
98) const
99{
100 const Foam::face& f = mesh.faces()[face()];
101
102 label faceBasePtI = mesh.tetBasePtIs()[face()];
103
104 if (faceBasePtI < 0)
105 {
106 faceBasePtI = 0;
107
108 if (warn && nWarnings_ < maxNWarnings)
109 {
110 ++nWarnings_;
112 << "No base point for face " << face() << ", " << f
113 << ", produces a valid tet decomposition." << endl;
114 if (nWarnings_ == maxNWarnings)
115 {
116 Warning
117 << "Suppressing any further warnings." << endl;
118 }
119 }
120 }
121
122 label facePtI = (tetPti_ + faceBasePtI) % f.size();
123 label faceOtherPtI = f.fcIndex(facePtI);
124
125 if (mesh.faceOwner()[face()] != cell())
126 {
127 std::swap(facePtI, faceOtherPtI);
128 }
129
130 return triFace(faceBasePtI, facePtI, faceOtherPtI);
131}
132
133
135{
136 const pointField& meshPoints = mesh.points();
137 const triFace tri = faceTriIs(mesh);
138
139 return tetPointRef
140 (
141 mesh.cellCentres()[cell()],
142 meshPoints[tri[0]],
143 meshPoints[tri[1]],
144 meshPoints[tri[2]]
145 );
146}
147
148
150{
151 const pointField& meshPoints = mesh.points();
152 const triFace tri = faceTriIs(mesh);
153
154 return triPointRef
155 (
156 meshPoints[tri[0]],
157 meshPoints[tri[1]],
158 meshPoints[tri[2]]
159 );
160}
161
162
164(
165 const polyMesh& mesh
166) const
167{
168 const pointField& meshOldPoints = mesh.oldPoints();
169 const triFace tri = faceTriIs(mesh);
170
171 return triPointRef
172 (
173 meshOldPoints[tri[0]],
174 meshOldPoints[tri[1]],
175 meshOldPoints[tri[2]]
176 );
177}
178
179
181(
182 const tetIndices& a,
183 const tetIndices& b
184) noexcept
185{
186 label diff;
187 return
188 (
189 ((diff = (a.cell() - b.cell())) != 0) ? diff
190 : ((diff = (a.face() - b.face())) != 0) ? diff
191 : (a.tetPt() - b.tetPt())
192 );
193}
194
195
196// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
197
198inline bool Foam::operator==(const tetIndices& a, const tetIndices& b) noexcept
199{
200 // Possibly slightly faster version than compare
201 return
202 (
203 a.cell() == b.cell()
204 && a.face() == b.face()
205 && a.tetPt() == b.tetPt()
206 );
207}
208
209
210inline bool Foam::operator!=(const tetIndices& a, const tetIndices& b) noexcept
211{
212 return !(a == b);
213}
214
215
216// ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label fcIndex(const label i) const noexcept
Definition: UListI.H:60
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:906
virtual const pointField & oldPoints() const
Return old points (mesh motion)
Definition: polyMesh.C:1133
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
const vectorField & cellCentres() const
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:84
triFace triIs(const polyMesh &mesh, const bool warn=true) const
Definition: tetIndicesI.H:95
static label compare(const tetIndices &a, const tetIndices &b) noexcept
Definition: tetIndicesI.H:181
triFace faceTriIs(const polyMesh &mesh, const bool warn=true) const
Return the indices corresponding to the tri on the face for.
Definition: tetIndicesI.H:55
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet.
Definition: tetIndicesI.H:134
triPointRef faceTri(const polyMesh &mesh) const
Definition: tetIndicesI.H:149
constexpr tetIndices() noexcept
Default construct, with invalid labels (-1)
Definition: tetIndicesI.H:31
triPointRef oldFaceTri(const polyMesh &mesh) const
Definition: tetIndicesI.H:164
A tetrahedron primitive.
Definition: tetrahedron.H:87
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:72
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:80
dynamicFvMesh & mesh
#define WarningInFunction
Report a warning using Foam::Warning.
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetPointRef.H:47
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangle.H:71
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
const direction noexcept
Definition: Scalar.H:223
messageStream Warning
labelList f(nPoints)
volScalarField & b
Definition: createFields.H:27