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