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 | Copyright (C) 2011-2016 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 
32 // Null constructor
34 :
35  index_(-3),
36  n_(Zero)
37 {}
38 
39 
40 // Construct from components
42 (
43  const label index,
44  const vector& n
45 )
46 :
47  index_(index),
48  n_(n)
49 {}
50 
51 
52 // Construct as copy
54 :
55  index_(w2.index()),
56  n_(w2.n())
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
62 template<class TrackingData>
63 inline bool Foam::directionInfo::valid(TrackingData& td) const
64 {
65  return index_ != -3;
66 }
67 
68 
69 // No geometric data so never any problem on cyclics
70 template<class TrackingData>
72 (
73  const polyMesh&,
74  const directionInfo& w2,
75  const scalar tol,
76  TrackingData& td
77 )
78  const
79 {
80  return true;
81 }
82 
83 
84 // index_ is already offset in face
85 template<class TrackingData>
87 (
88  const polyMesh&,
89  const polyPatch& patch,
90  const label patchFacei,
91  const point& faceCentre,
92  TrackingData& td
93 )
94 {}
95 
96 
97 // index_ is offset in face on other side. So reverse it here.
98 // (Note: f[0] on other domain is connected to f[0] in this domain,
99 // f[1] ,, f[size-1] ,,
100 // etc.)
101 template<class TrackingData>
103 (
104  const polyMesh&,
105  const polyPatch& patch,
106  const label patchFacei,
107  const point& faceCentre,
108  TrackingData& td
109 )
110 {
111  if (index_ >= 0)
112  {
113  const face& f = patch[patchFacei];
114 
115  index_ = (f.size() - index_) % f.size();
116  }
117 }
118 
119 
120 // No geometric data.
121 template<class TrackingData>
123 (
124  const polyMesh&,
125  const tensor& rotTensor,
126  TrackingData& td
127 )
128 {}
129 
130 
131 // Update this cell with neighbouring face information
132 template<class TrackingData>
134 (
135  const polyMesh& mesh,
136  const label thisCelli,
137  const label neighbourFacei,
138  const directionInfo& neighbourInfo,
139  const scalar, // tol
140  TrackingData& td
141 )
142 {
143  if (index_ >= -2)
144  {
145  // Already determined.
146  return false;
147  }
148 
149  if (hexMatcher().isA(mesh, thisCelli))
150  {
151  const face& f = mesh.faces()[neighbourFacei];
152 
153  if (neighbourInfo.index() == -2)
154  {
155  // Geometric information from neighbour
156  index_ = -2;
157  }
158  else if (neighbourInfo.index() == -1)
159  {
160  // Cut tangential to face. Take any edge connected to face
161  // but not used in face.
162 
163  // Get first edge on face.
164  label edgeI = mesh.faceEdges()[neighbourFacei][0];
165 
166  const edge& e = mesh.edges()[edgeI];
167 
168  // Find face connected to face through edgeI and on same cell.
169  label facei =
171  (
172  mesh,
173  thisCelli,
174  neighbourFacei,
175  edgeI
176  );
177 
178  // Find edge on facei which is connected to e.start() but not edgeI.
179  index_ =
181  (
182  mesh,
183  mesh.faceEdges()[facei],
184  edgeI,
185  e.start()
186  );
187  }
188  else
189  {
190  // Index is a vertex on the face. Convert to mesh edge.
191 
192  // Get mesh edge between f[index_] and f[index_+1]
193  label v0 = f[neighbourInfo.index()];
194  label v1 = f[(neighbourInfo.index() + 1) % f.size()];
195 
196  index_ = findEdge(mesh, mesh.faceEdges()[neighbourFacei], v0, v1);
197  }
198  }
199  else
200  {
201  // Not a hex so mark this as geometric.
202  index_ = -2;
203  }
204 
205 
206  n_ = neighbourInfo.n();
207 
208  return true;
209 }
210 
211 
212 // Update this face with neighbouring cell information
213 template<class TrackingData>
215 (
216  const polyMesh& mesh,
217  const label thisFacei,
218  const label neighbourCelli,
219  const directionInfo& neighbourInfo,
220  const scalar, // tol
221  TrackingData& td
222 )
223 {
224  // Handle special cases first
225 
226  if (index_ >= -2)
227  {
228  // Already determined
229  return false;
230  }
231 
232  // Handle normal cases where topological or geometrical info comes from
233  // neighbouring cell
234 
235  if (neighbourInfo.index() >= 0)
236  {
237  // Neighbour has topological direction (and hence is hex). Find cut
238  // edge on face.
239  index_ =
241  (
242  mesh,
243  neighbourCelli,
244  thisFacei,
245  neighbourInfo.index()
246  );
247  }
248  else
249  {
250  // Neighbour has geometric information. Use.
251  index_ = -2;
252  }
253 
254 
255  n_ = neighbourInfo.n();
256 
257  return true;
258 }
259 
260 
261 // Merge this with information on same face
262 template<class TrackingData>
264 (
265  const polyMesh& mesh,
266  const label, // thisFacei
267  const directionInfo& neighbourInfo,
268  const scalar, // tol
269  TrackingData& td
270 )
271 {
272  if (index_ >= -2)
273  {
274  // Already visited.
275  return false;
276  }
277  else
278  {
279  index_ = neighbourInfo.index();
280 
281  n_ = neighbourInfo.n();
282 
283  return true;
284  }
285 }
286 
287 
288 template<class TrackingData>
289 inline bool Foam::directionInfo::equal
290 (
291  const directionInfo& rhs,
292  TrackingData& td
293 ) const
294 {
295  return operator==(rhs);
296 }
297 
298 
299 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
300 
302  const
303 {
304  return index() == rhs.index() && n() == rhs.n();
305 }
306 
307 
309  const
310 {
311  return !(*this == rhs);
312 }
313 
314 
315 // ************************************************************************* //
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
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
const double e
Elementary charge.
Definition: doubleFloat.H:78
bool equal(const directionInfo &, TrackingData &td) const
Same (like operator==)
const vector & n() const
label index() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
A cellMatcher for hex cells.
Definition: hexMatcher.H:51
void leaveDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Convert any absolute coordinates into relative to (patch)face.
bool updateCell(const polyMesh &, const label thisCelli, const label neighbourFacei, const directionInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.
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
void enterDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Reverse of leaveDomain.
bool operator!=(const directionInfo &) const
bool operator==(const directionInfo &) const
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:91
label start() const
Return start vertex label.
Definition: edgeI.H:81
labelList f(nPoints)
Holds direction in which to split cell (in fact a local coordinate axes). Information is a label and ...
Definition: directionInfo.H:84
directionInfo()
Construct null.
void transform(const polyMesh &, const tensor &, TrackingData &td)
Apply rotation matrix to any coordinates.
bool updateFace(const polyMesh &, const label thisFacei, const label neighbourCelli, const directionInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring cell.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
const labelListList & faceEdges() const
bool sameGeometry(const polyMesh &, const directionInfo &, const scalar, TrackingData &td) const
Check for identical geometrical data. Used for cyclics checking.
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