directionInfoI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "polyMesh.H"
27 #include "meshTools.H"
28 #include "hexMatcher.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
33 :
34  index_(-3),
35  n_(Zero)
36 {}
37 
38 
40 (
41  const label index,
42  const vector& n
43 )
44 :
45  index_(index),
46  n_(n)
47 {}
48 
49 
50 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
51 
52 template<class TrackingData>
53 inline bool Foam::directionInfo::valid(TrackingData& td) const
54 {
55  return index_ != -3;
56 }
57 
58 
59 template<class TrackingData>
61 (
62  const polyMesh&,
63  const directionInfo& w2,
64  const scalar tol,
65  TrackingData& td
66 )
67  const
68 {
69  return true;
70 }
71 
72 
73 template<class TrackingData>
75 (
76  const polyPatch& patch,
77  const label patchFacei,
78  const transformer& transform,
79  TrackingData& td
80 )
81 {
82  if (index_ >= 0)
83  {
84  const face& f = patch[patchFacei];
85  index_ = (f.size() - index_) % f.size();
86  }
87 }
88 
89 
90 template<class TrackingData>
92 (
93  const polyMesh& mesh,
94  const label thisCelli,
95  const label neighbourFacei,
96  const directionInfo& neighbourInfo,
97  const scalar, // tol
98  TrackingData& td
99 )
100 {
101  if (index_ >= -2)
102  {
103  // Already determined.
104  return false;
105  }
106 
107  if (hexMatcher().isA(mesh, thisCelli))
108  {
109  const face& f = mesh.faces()[neighbourFacei];
110 
111  if (neighbourInfo.index() == -2)
112  {
113  // Geometric information from neighbour
114  index_ = -2;
115  }
116  else if (neighbourInfo.index() == -1)
117  {
118  // Cut tangential to face. Take any edge connected to face
119  // but not used in face.
120 
121  // Get first edge on face.
122  label edgeI = mesh.faceEdges()[neighbourFacei][0];
123 
124  const edge& e = mesh.edges()[edgeI];
125 
126  // Find face connected to face through edgeI and on same cell.
127  label facei =
129  (
130  mesh,
131  thisCelli,
132  neighbourFacei,
133  edgeI
134  );
135 
136  // Find edge on facei which is connected to e.start() but not edgeI.
137  index_ =
139  (
140  mesh,
141  mesh.faceEdges()[facei],
142  edgeI,
143  e.start()
144  );
145  }
146  else
147  {
148  // Index is a vertex on the face. Convert to mesh edge.
149 
150  // Get mesh edge between f[index_] and f[index_+1]
151  label v0 = f[neighbourInfo.index()];
152  label v1 = f[(neighbourInfo.index() + 1) % f.size()];
153 
154  index_ = findEdge(mesh, mesh.faceEdges()[neighbourFacei], v0, v1);
155  }
156  }
157  else
158  {
159  // Not a hex so mark this as geometric.
160  index_ = -2;
161  }
162 
163 
164  n_ = neighbourInfo.n();
165 
166  return true;
167 }
168 
169 
170 template<class TrackingData>
172 (
173  const polyMesh& mesh,
174  const label thisFacei,
175  const label neighbourCelli,
176  const directionInfo& neighbourInfo,
177  const scalar, // tol
178  TrackingData& td
179 )
180 {
181  // Handle special cases first
182 
183  if (index_ >= -2)
184  {
185  // Already determined
186  return false;
187  }
188 
189  // Handle normal cases where topological or geometrical info comes from
190  // neighbouring cell
191 
192  if (neighbourInfo.index() >= 0)
193  {
194  // Neighbour has topological direction (and hence is hex). Find cut
195  // edge on face.
196  index_ =
198  (
199  mesh,
200  neighbourCelli,
201  thisFacei,
202  neighbourInfo.index()
203  );
204  }
205  else
206  {
207  // Neighbour has geometric information. Use.
208  index_ = -2;
209  }
210 
211 
212  n_ = neighbourInfo.n();
213 
214  return true;
215 }
216 
217 
218 template<class TrackingData>
220 (
221  const polyMesh& mesh,
222  const label, // thisFacei
223  const directionInfo& neighbourInfo,
224  const scalar, // tol
225  TrackingData& td
226 )
227 {
228  if (index_ >= -2)
229  {
230  // Already visited.
231  return false;
232  }
233  else
234  {
235  index_ = neighbourInfo.index();
236 
237  n_ = neighbourInfo.n();
238 
239  return true;
240  }
241 }
242 
243 
244 template<class TrackingData>
245 inline bool Foam::directionInfo::equal
246 (
247  const directionInfo& rhs,
248  TrackingData& td
249 ) const
250 {
251  return operator==(rhs);
252 }
253 
254 
255 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
256 
257 inline bool Foam::directionInfo::operator==
258 (
259  const Foam::directionInfo& rhs
260 ) const
261 {
262  return index() == rhs.index() && n() == rhs.n();
263 }
264 
265 
266 inline bool Foam::directionInfo::operator!=
267 (
268  const Foam::directionInfo& rhs
269 ) const
270 {
271  return !(*this == rhs);
272 }
273 
274 
275 // ************************************************************************* //
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
bool sameGeometry(const polyMesh &, const directionInfo &, const scalar, TrackingData &td) const
Check for identical geometrical data. Used for cyclics checking.
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
label otherEdge(const primitiveMesh &, const labelList &edgeLabels, const label edgeI, const label vertI)
Return label of other edge (out of candidates edgeLabels)
Definition: meshTools.C:501
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
bool equal(const directionInfo &, TrackingData &td) const
Same (like operator==)
const labelListList & faceEdges() const
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
bool operator==(const directionInfo &) const
A cellMatcher for hex cells.
Definition: hexMatcher.H:51
bool updateCell(const polyMesh &, const label thisCelli, const label neighbourFacei, const directionInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.
void transform(const polyPatch &patch, const label patchFacei, const transformer &transform, TrackingData &td)
Transform across an interface.
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
static label edgeToFaceIndex(const primitiveMesh &mesh, const label celli, const label facei, const label edgeI)
Given edge on hex cell find corresponding edge on face. Is either.
Definition: directionInfo.C:93
static const zero Zero
Definition: zero.H:97
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
label index() const
labelList f(nPoints)
Holds direction in which to split cell (in fact a local coordinate axes). Information is a label and ...
Definition: directionInfo.H:85
directionInfo()
Construct null.
const vector & n() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
bool updateFace(const polyMesh &, const label thisFacei, const label neighbourCelli, const directionInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring cell.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
label start() const
Return start vertex label.
Definition: edgeI.H:81
label otherFace(const primitiveMesh &, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:536