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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30 
31 inline Foam::label Foam::tetIndices::cell() const
32 {
33  return celli_;
34 }
35 
36 
37 inline Foam::label& Foam::tetIndices::cell()
38 {
39  return celli_;
40 }
41 
42 
43 inline Foam::label Foam::tetIndices::face() const
44 {
45  return facei_;
46 }
47 
48 
49 inline Foam::label& Foam::tetIndices::face()
50 {
51  return facei_;
52 }
53 
54 
55 inline Foam::label Foam::tetIndices::tetPt() const
56 {
57  return tetPti_;
58 }
59 
60 
61 inline Foam::label& Foam::tetIndices::tetPt()
62 {
63  return tetPti_;
64 }
65 
66 
68 (
69  const polyMesh& mesh,
70  const bool warn
71 ) const
72 {
73  const Foam::face& f = mesh.faces()[face()];
74 
75  label faceBasePtI = mesh.tetBasePtIs()[face()];
76 
77  if (faceBasePtI < 0)
78  {
79  faceBasePtI = 0;
80 
81  if (warn)
82  {
83  if (nWarnings < maxNWarnings)
84  {
86  << "No base point for face " << face() << ", " << f
87  << ", produces a valid tet decomposition." << endl;
88  ++nWarnings;
89  }
90  if (nWarnings == maxNWarnings)
91  {
92  Warning
93  << "Suppressing any further warnings." << endl;
94  ++nWarnings;
95  }
96  }
97  }
98 
99  label facePtI = (tetPti_ + faceBasePtI) % f.size();
100  label faceOtherPtI = f.fcIndex(facePtI);
101 
102  if (mesh.faceOwner()[face()] != cell())
103  {
104  std::swap(facePtI, faceOtherPtI);
105  }
106 
107  return triFace(f[faceBasePtI], f[facePtI], f[faceOtherPtI]);
108 }
109 
110 
112 (
113  const polyMesh& mesh,
114  const bool warn
115 ) const
116 {
117  const Foam::face& f = mesh.faces()[face()];
118 
119  label faceBasePtI = mesh.tetBasePtIs()[face()];
120 
121  if (faceBasePtI < 0)
122  {
123  faceBasePtI = 0;
124 
125  if (warn)
126  {
127  if (nWarnings < maxNWarnings)
128  {
130  << "No base point for face " << face() << ", " << f
131  << ", produces a valid tet decomposition." << endl;
132  ++nWarnings;
133  }
134  if (nWarnings == maxNWarnings)
135  {
136  Warning
137  << "Suppressing any further warnings." << endl;
138  ++nWarnings;
139  }
140  }
141  }
142 
143  label facePtI = (tetPti_ + faceBasePtI) % f.size();
144  label faceOtherPtI = f.fcIndex(facePtI);
145 
146  if (mesh.faceOwner()[face()] != cell())
147  {
148  std::swap(facePtI, faceOtherPtI);
149  }
150 
151  return triFace(faceBasePtI, facePtI, faceOtherPtI);
152 }
153 
154 
156 {
157  const pointField& meshPoints = mesh.points();
158  const triFace tri = faceTriIs(mesh);
159 
160  return tetPointRef
161  (
162  mesh.cellCentres()[cell()],
163  meshPoints[tri[0]],
164  meshPoints[tri[1]],
165  meshPoints[tri[2]]
166  );
167 }
168 
169 
171 {
172  const pointField& meshPoints = mesh.points();
173  const triFace tri = faceTriIs(mesh);
174 
175  return triPointRef
176  (
177  meshPoints[tri[0]],
178  meshPoints[tri[1]],
179  meshPoints[tri[2]]
180  );
181 }
182 
183 
185 (
186  const polyMesh& mesh
187 ) const
188 {
189  const pointField& meshOldPoints = mesh.oldPoints();
190  const triFace tri = faceTriIs(mesh);
191 
192  return triPointRef
193  (
194  meshOldPoints[tri[0]],
195  meshOldPoints[tri[1]],
196  meshOldPoints[tri[2]]
197  );
198 }
199 
200 
201 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
202 
203 inline bool Foam::tetIndices::operator==(const Foam::tetIndices& rhs) const
204 {
205  return
206  cell() == rhs.cell()
207  && face() == rhs.face()
208  && tetPt() == rhs.tetPt();
209 }
210 
211 
212 inline bool Foam::tetIndices::operator!=(const Foam::tetIndices& rhs) const
213 {
214  return !(*this == rhs);
215 }
216 
217 
218 // ************************************************************************* //
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::tetIndices::tet
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet.
Definition: tetIndicesI.H:155
Foam::tetIndices::triIs
triFace triIs(const polyMesh &mesh, const bool warn=true) const
Return the local indices corresponding to the tri on the face.
Definition: tetIndicesI.H:112
Foam::tetIndices::operator==
bool operator==(const tetIndices &) const
Definition: tetIndicesI.H:203
Foam::tetIndices::cell
label cell() const
Return the cell.
Definition: tetIndicesI.H:31
Foam::Warning
messageStream Warning
Foam::tetIndices::face
label face() const
Return the face.
Definition: tetIndicesI.H:43
Foam::tetIndices::oldFaceTri
triPointRef oldFaceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
Definition: tetIndicesI.H:185
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::tetPointRef
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetPointRef.H:47
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:59
Foam::Field< vector >
Foam::tetIndices::faceTri
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
Definition: tetIndicesI.H:170
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::polyMesh::oldPoints
virtual const pointField & oldPoints() const
Return old points (mesh motion)
Definition: polyMesh.C:1119
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::tetIndices::faceTriIs
triFace faceTriIs(const polyMesh &mesh, const bool warn=true) const
Return the indices corresponding to the tri on the face for.
Definition: tetIndicesI.H:68
Foam::tetIndices::tetPt
label tetPt() const
Return the characterising tetPtI.
Definition: tetIndicesI.H:55
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:83
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:69
f
labelList f(nPoints)
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::tetIndices::operator!=
bool operator!=(const tetIndices &) const
Definition: tetIndicesI.H:212
triFace
face triFace(3)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::triPointRef
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangle.H:71
Foam::tetrahedron
A tetrahedron primitive.
Definition: tetrahedron.H:63
Foam::polyMesh::tetBasePtIs
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:892