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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 inline Foam::label Foam::tetIndices::cell() const
31 {
32  return celli_;
33 }
34 
35 
36 inline Foam::label& Foam::tetIndices::cell()
37 {
38  return celli_;
39 }
40 
41 
42 inline Foam::label Foam::tetIndices::face() const
43 {
44  return facei_;
45 }
46 
47 
48 inline Foam::label& Foam::tetIndices::face()
49 {
50  return facei_;
51 }
52 
53 
54 inline Foam::label Foam::tetIndices::tetPt() const
55 {
56  return tetPti_;
57 }
58 
59 
60 inline Foam::label& Foam::tetIndices::tetPt()
61 {
62  return tetPti_;
63 }
64 
65 
67 (
68  const polyMesh& mesh,
69  const bool warn
70 ) const
71 {
72  const Foam::face& f = mesh.faces()[face()];
73 
74  label faceBasePtI = mesh.tetBasePtIs()[face()];
75 
76  if (faceBasePtI < 0)
77  {
78  faceBasePtI = 0;
79 
80  if (warn)
81  {
82  if (nWarnings < maxNWarnings)
83  {
85  << "No base point for face " << face() << ", " << f
86  << ", produces a valid tet decomposition." << endl;
87  ++nWarnings;
88  }
89  if (nWarnings == maxNWarnings)
90  {
91  Warning
92  << "Suppressing any further warnings." << endl;
93  ++nWarnings;
94  }
95  }
96  }
97 
98  label facePtI = (tetPti_ + faceBasePtI) % f.size();
99  label faceOtherPtI = f.fcIndex(facePtI);
100 
101  if (mesh.faceOwner()[face()] != cell())
102  {
103  Swap(facePtI, faceOtherPtI);
104  }
105 
106  return triFace(f[faceBasePtI], f[facePtI], f[faceOtherPtI]);
107 }
108 
109 
111 (
112  const polyMesh& mesh,
113  const bool warn
114 ) const
115 {
116  const Foam::face& f = mesh.faces()[face()];
117 
118  label faceBasePtI = mesh.tetBasePtIs()[face()];
119 
120  if (faceBasePtI < 0)
121  {
122  faceBasePtI = 0;
123 
124  if (warn)
125  {
126  if (nWarnings < maxNWarnings)
127  {
129  << "No base point for face " << face() << ", " << f
130  << ", produces a valid tet decomposition." << endl;
131  ++nWarnings;
132  }
133  if (nWarnings == maxNWarnings)
134  {
135  Warning
136  << "Suppressing any further warnings." << endl;
137  ++nWarnings;
138  }
139  }
140  }
141 
142  label facePtI = (tetPti_ + faceBasePtI) % f.size();
143  label faceOtherPtI = f.fcIndex(facePtI);
144 
145  if (mesh.faceOwner()[face()] != cell())
146  {
147  Swap(facePtI, faceOtherPtI);
148  }
149 
150  return triFace(faceBasePtI, facePtI, faceOtherPtI);
151 }
152 
153 
155 {
156  const pointField& meshPoints = mesh.points();
157  const triFace tri = faceTriIs(mesh);
158 
159  return tetPointRef
160  (
161  mesh.cellCentres()[cell()],
162  meshPoints[tri[0]],
163  meshPoints[tri[1]],
164  meshPoints[tri[2]]
165  );
166 }
167 
168 
170 {
171  const pointField& meshPoints = mesh.points();
172  const triFace tri = faceTriIs(mesh);
173 
174  return triPointRef
175  (
176  meshPoints[tri[0]],
177  meshPoints[tri[1]],
178  meshPoints[tri[2]]
179  );
180 }
181 
182 
184 (
185  const polyMesh& mesh
186 ) const
187 {
188  const pointField& meshOldPoints = mesh.oldPoints();
189  const triFace tri = faceTriIs(mesh);
190 
191  return triPointRef
192  (
193  meshOldPoints[tri[0]],
194  meshOldPoints[tri[1]],
195  meshOldPoints[tri[2]]
196  );
197 }
198 
199 
200 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
201 
202 inline bool Foam::tetIndices::operator==(const Foam::tetIndices& rhs) const
203 {
204  return
205  cell() == rhs.cell()
206  && face() == rhs.face()
207  && tetPt() == rhs.tetPt();
208 }
209 
210 
211 inline bool Foam::tetIndices::operator!=(const Foam::tetIndices& rhs) const
212 {
213  return !(*this == rhs);
214 }
215 
216 
217 // ************************************************************************* //
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:154
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:111
Foam::tetIndices::operator==
bool operator==(const tetIndices &) const
Definition: tetIndicesI.H:202
Foam::tetIndices::cell
label cell() const
Return the cell.
Definition: tetIndicesI.H:30
Foam::Warning
messageStream Warning
Foam::tetIndices::face
label face() const
Return the face.
Definition: tetIndicesI.H:42
Foam::tetIndices::oldFaceTri
triPointRef oldFaceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
Definition: tetIndicesI.H:184
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::Swap
void Swap(DynamicList< T, SizeMin1 > &a, DynamicList< T, SizeMin2 > &b)
Definition: DynamicListI.H:913
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:61
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:169
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:67
Foam::tetIndices::tetPt
label tetPt() const
Return the characterising tetPtI.
Definition: tetIndicesI.H:54
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:211
triFace
face triFace(3)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
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:73
Foam::tetrahedron
A tetrahedron primitive.
Definition: tetrahedron.H:65
Foam::polyMesh::tetBasePtIs
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:892