directionInfoI.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-2016 OpenFOAM Foundation
9  Copyright (C) 2020 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 #include "polyMesh.H"
30 #include "meshTools.H"
31 #include "hexMatcher.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 :
37  index_(-3),
38  n_(Zero)
39 {}
40 
41 
43 (
44  const label index,
45  const vector& n
46 )
47 :
48  index_(index),
49  n_(n)
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
54 
55 template<class TrackingData>
56 inline bool Foam::directionInfo::valid(TrackingData& td) const
57 {
58  return index_ != -3;
59 }
60 
61 
62 // No geometric data so never any problem on cyclics
63 template<class TrackingData>
65 (
66  const polyMesh&,
67  const directionInfo& w2,
68  const scalar tol,
69  TrackingData& td
70 ) const
71 {
72  return true;
73 }
74 
75 
76 // index_ is already offset in face
77 template<class TrackingData>
79 (
80  const polyMesh&,
81  const polyPatch& patch,
82  const label patchFacei,
83  const point& faceCentre,
84  TrackingData& td
85 )
86 {}
87 
88 
89 // index_ is offset in face on other side. So reverse it here.
90 // (Note: f[0] on other domain is connected to f[0] in this domain,
91 // f[1] ,, f[size-1] ,,
92 // etc.)
93 template<class TrackingData>
95 (
96  const polyMesh&,
97  const polyPatch& patch,
98  const label patchFacei,
99  const point& faceCentre,
100  TrackingData& td
101 )
102 {
103  if (index_ >= 0)
104  {
105  const face& f = patch[patchFacei];
106 
107  index_ = (f.size() - index_) % f.size();
108  }
109 }
110 
111 
112 // No geometric data.
113 template<class TrackingData>
115 (
116  const polyMesh&,
117  const tensor& rotTensor,
118  TrackingData& td
119 )
120 {}
121 
122 
123 // Update this cell with neighbouring face information
124 template<class TrackingData>
126 (
127  const polyMesh& mesh,
128  const label thisCelli,
129  const label neighbourFacei,
130  const directionInfo& neighbourInfo,
131  const scalar, // tol
132  TrackingData& td
133 )
134 {
135  if (index_ >= -2)
136  {
137  // Already determined.
138  return false;
139  }
140 
141  if (hexMatcher::test(mesh, thisCelli))
142  {
143  const face& f = mesh.faces()[neighbourFacei];
144 
145  if (neighbourInfo.index() == -2)
146  {
147  // Geometric information from neighbour
148  index_ = -2;
149  }
150  else if (neighbourInfo.index() == -1)
151  {
152  // Cut tangential to face. Take any edge connected to face
153  // but not used in face.
154 
155  // Get first edge on face.
156  label edgeI = mesh.faceEdges()[neighbourFacei][0];
157 
158  const edge& e = mesh.edges()[edgeI];
159 
160  // Find face connected to face through edgeI and on same cell.
161  label facei =
163  (
164  mesh,
165  thisCelli,
166  neighbourFacei,
167  edgeI
168  );
169 
170  // Find edge on facei which is connected to e.start() but not edgeI.
171  index_ =
173  (
174  mesh,
175  mesh.faceEdges()[facei],
176  edgeI,
177  e.start()
178  );
179  }
180  else
181  {
182  // Index is a vertex on the face. Convert to mesh edge.
183 
184  // Get mesh edge between f[index_] and f[index_+1]
185  label v0 = f[neighbourInfo.index()];
186  label v1 = f[(neighbourInfo.index() + 1) % f.size()];
187 
188  index_ = findEdge(mesh, mesh.faceEdges()[neighbourFacei], v0, v1);
189  }
190  }
191  else
192  {
193  // Not a hex so mark this as geometric.
194  index_ = -2;
195  }
196 
197 
198  n_ = neighbourInfo.n();
199 
200  return true;
201 }
202 
203 
204 // Update this face with neighbouring cell information
205 template<class TrackingData>
207 (
208  const polyMesh& mesh,
209  const label thisFacei,
210  const label neighbourCelli,
211  const directionInfo& neighbourInfo,
212  const scalar, // tol
213  TrackingData& td
214 )
215 {
216  // Handle special cases first
217 
218  if (index_ >= -2)
219  {
220  // Already determined
221  return false;
222  }
223 
224  // Handle normal cases where topological or geometrical info comes from
225  // neighbouring cell
226 
227  if (neighbourInfo.index() >= 0)
228  {
229  // Neighbour has topological direction (and hence is hex). Find cut
230  // edge on face.
231  index_ =
232  edgeToFaceIndex
233  (
234  mesh,
235  neighbourCelli,
236  thisFacei,
237  neighbourInfo.index()
238  );
239  }
240  else
241  {
242  // Neighbour has geometric information. Use.
243  index_ = -2;
244  }
245 
246  n_ = neighbourInfo.n();
247 
248  return true;
249 }
250 
251 
252 // Merge this with information on same face
253 template<class TrackingData>
255 (
256  const polyMesh& mesh,
257  const label, // thisFacei
258  const directionInfo& neighbourInfo,
259  const scalar, // tol
260  TrackingData& td
261 )
262 {
263  if (index_ >= -2)
264  {
265  // Already visited.
266  return false;
267  }
268  else
269  {
270  index_ = neighbourInfo.index();
271 
272  n_ = neighbourInfo.n();
273 
274  return true;
275  }
276 }
277 
278 
279 template<class TrackingData>
280 inline bool Foam::directionInfo::equal
281 (
282  const directionInfo& rhs,
283  TrackingData& td
284 ) const
285 {
286  return operator==(rhs);
287 }
288 
289 
290 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
291 
292 inline bool Foam::directionInfo::operator==
293 (
294  const directionInfo& rhs
295 ) const
296 {
297  return index_ == rhs.index_ && n_ == rhs.n_;
298 }
299 
300 
301 inline bool Foam::directionInfo::operator!=
302 (
303  const directionInfo& rhs
304 ) const
305 {
306  return !(*this == rhs);
307 }
308 
309 
310 // ************************************************************************* //
Foam::Tensor< scalar >
meshTools.H
Foam::directionInfo::index
label index() const
Definition: directionInfo.H:146
Foam::meshTools::otherFace
label otherFace(const primitiveMesh &mesh, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:555
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::directionInfo::enterDomain
void enterDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Reverse of leaveDomain.
Definition: directionInfoI.H:95
hexMatcher.H
Foam::directionInfo::equal
bool equal(const directionInfo &, TrackingData &td) const
Test for equality, with TrackingData.
Definition: directionInfoI.H:281
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::hexMatcher::test
static bool test(const UList< face > &faces)
Definition: hexMatcher.C:87
Foam::directionInfo::updateFace
bool updateFace(const polyMesh &, const label thisFacei, const label neighbourCelli, const directionInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring cell.
Definition: directionInfoI.H:207
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::directionInfo::valid
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
Definition: directionInfoI.H:56
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::directionInfo::transform
void transform(const polyMesh &, const tensor &, TrackingData &td)
Apply rotation matrix to any coordinates.
Definition: directionInfoI.H:115
Foam::meshTools::findEdge
label findEdge(const edgeList &edges, const labelList &candidates, const label v0, const label v1)
Return edge among candidates that uses the two vertices.
Definition: meshTools.C:359
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::primitiveMesh::faceEdges
const labelListList & faceEdges() const
Definition: primitiveMeshEdges.C:528
Foam::directionInfo
Holds direction in which to split cell (in fact a local coordinate axes). Information is a label and ...
Definition: directionInfo.H:83
w2
#define w2
Definition: blockCreate.C:35
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::meshTools::otherEdge
label otherEdge(const primitiveMesh &mesh, const labelList &edgeLabels, const label thisEdgeI, const label thisVertI)
Return label of other edge (out of candidates edgeLabels)
Definition: meshTools.C:521
Foam::directionInfo::n
const vector & n() const
Definition: directionInfo.H:151
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::directionInfo::leaveDomain
void leaveDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Convert any absolute coordinates into relative to (patch)face.
Definition: directionInfoI.H:79
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::directionInfo::directionInfo
directionInfo()
Default construct, index=-1, vector::zero.
Definition: directionInfoI.H:35
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::directionInfo::sameGeometry
bool sameGeometry(const polyMesh &, const directionInfo &, const scalar, TrackingData &td) const
Check for identical geometrical data (eg, cyclics checking)
Definition: directionInfoI.H:65
Foam::directionInfo::updateCell
bool updateCell(const polyMesh &, const label thisCelli, const label neighbourFacei, const directionInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.
Definition: directionInfoI.H:126