BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / dynamicMesh / meshCut / directions / directionInfo / directionInfoI.H
blob2f0fdd004f7c8f6c07fdf45be636d1fd5110ae55
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "polyMesh.H"
27 #include "meshTools.H"
28 #include "hexMatcher.H"
30 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
32 // Null constructor
33 inline Foam::directionInfo::directionInfo()
35     index_(-3),
36     n_(vector::zero)
40 // Construct from components
41 inline Foam::directionInfo::directionInfo
43     const label index,
44     const vector& n
47     index_(index),
48     n_(n)
52 // Construct as copy
53 inline Foam::directionInfo::directionInfo(const directionInfo& w2)
55     index_(w2.index()),
56     n_(w2.n())
60 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
62 template<class TrackingData>
63 inline bool Foam::directionInfo::valid(TrackingData& td) const
65     return index_ != -3;
69 // No geometric data so never any problem on cyclics
70 template<class TrackingData>
71 inline bool Foam::directionInfo::sameGeometry
73     const polyMesh&,
74     const directionInfo& w2,
75     const scalar tol,
76     TrackingData& td
78  const
80     return true;
84 // index_ is already offset in face
85 template<class TrackingData>
86 inline void Foam::directionInfo::leaveDomain
88     const polyMesh&,
89     const polyPatch& patch,
90     const label patchFaceI,
91     const point& faceCentre,
92     TrackingData& td
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>
102 inline void Foam::directionInfo::enterDomain
104     const polyMesh&,
105     const polyPatch& patch,
106     const label patchFaceI,
107     const point& faceCentre,
108     TrackingData& td
111     if (index_ >= 0)
112     {
113         const face& f = patch[patchFaceI];
115         index_ = (f.size() - index_) % f.size();
116     }
120 // No geometric data.
121 template<class TrackingData>
122 inline void Foam::directionInfo::transform
124     const polyMesh&,
125     const tensor& rotTensor,
126     TrackingData& td
131 // Update this cell with neighbouring face information
132 template<class TrackingData>
133 inline bool Foam::directionInfo::updateCell
135     const polyMesh& mesh,
136     const label thisCellI,
137     const label neighbourFaceI,
138     const directionInfo& neighbourInfo,
139     const scalar,       // tol
140     TrackingData& td
143     if (index_ >= -2)
144     {
145         // Already determined.
146         return false;
147     }
149     if (hexMatcher().isA(mesh, thisCellI))
150     {
151         const face& f = mesh.faces()[neighbourFaceI];
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.
163             // Get first edge on face.
164             label edgeI = mesh.faceEdges()[neighbourFaceI][0];
166             const edge& e = mesh.edges()[edgeI];
168             // Find face connected to face through edgeI and on same cell.
169             label faceI =
170                 meshTools::otherFace
171                 (
172                     mesh,
173                     thisCellI,
174                     neighbourFaceI,
175                     edgeI
176                 );
178             // Find edge on faceI which is connected to e.start() but not edgeI.
179             index_ =
180                 meshTools::otherEdge
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.
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()];
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     }
206     n_ = neighbourInfo.n();
208     return true;
212 // Update this face with neighbouring cell information
213 template<class TrackingData>
214 inline bool Foam::directionInfo::updateFace
216     const polyMesh& mesh,
217     const label thisFaceI,
218     const label neighbourCellI,
219     const directionInfo& neighbourInfo,
220     const scalar,   // tol
221     TrackingData& td
224     // Handle special cases first
226     if (index_ >= -2)
227     {
228         // Already determined
229         return false;
230     }
232     // Handle normal cases where topological or geometrical info comes from
233     // neighbouring cell
235     if (neighbourInfo.index() >= 0)
236     {
237         // Neighbour has topological direction (and hence is hex). Find cut
238         // edge on face.
239         index_ =
240             edgeToFaceIndex
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     }
255     n_ = neighbourInfo.n();
257     return true;
261 // Merge this with information on same face
262 template<class TrackingData>
263 inline bool Foam::directionInfo::updateFace
265     const polyMesh& mesh,
266     const label,    // thisFaceI
267     const directionInfo& neighbourInfo,
268     const scalar,   // tol
269     TrackingData& td
272     if (index_ >= -2)
273     {
274         // Already visited.
275         return false;
276     }
277     else
278     {
279         index_ = neighbourInfo.index();
281         n_ = neighbourInfo.n();
283         return true;
284     }
288 template <class TrackingData>
289 inline bool Foam::directionInfo::equal
291     const directionInfo& rhs,
292     TrackingData& td
293 ) const
295     return operator==(rhs);
299 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
301 inline bool Foam::directionInfo::operator==(const Foam::directionInfo& rhs)
302  const
304     return  index() == rhs.index() && n() == rhs.n();
308 inline bool Foam::directionInfo::operator!=(const Foam::directionInfo& rhs)
309  const
311     return !(*this == rhs);
315 // ************************************************************************* //