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