BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / dynamicMesh / meshCut / directions / directionInfo / directionInfo.H
blob88d2493302da230b66b2a176edf16bbd46b165a0
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 Class
25     Foam::directionInfo
27 Description
28     Holds direction in which to split cell (in fact a local coordinate axes).
29     Information is a label and a direction.
31     The direction is the normal
32     direction to cut in. The label's meaning depends on whether the info
33     is on a cell or on a face:
34         - in cell: edge that is being cut. (determines for hex how cut is)
35         - in face: local face point that is being cut or -1.
36             -# (-1)  : cut is tangential to plane
37             -# (>= 0): edge fp..fp+1 is cut
39             (has to be facepoint, not vertex since vertex not valid across
40              processors whereas f[0] should correspond to f[0] on other side)
42     The rule is that if the label is set (-1 or higher) it is used
43     (topological information only), otherwise the vector is used. This makes
44     sure that we use topological information as much as possible and so a
45     hex mesh is cut purely topologically. All other shapes are cut
46     geometrically.
48 SourceFiles
49     directionInfoI.H
50     directionInfo.C
52 \*---------------------------------------------------------------------------*/
54 #ifndef directionInfo_H
55 #define directionInfo_H
57 #include "point.H"
58 #include "labelList.H"
59 #include "tensor.H"
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 namespace Foam
65 class polyPatch;
66 class polyMesh;
67 class primitiveMesh;
68 class edge;
69 class face;
70 class polyMesh;
72 /*---------------------------------------------------------------------------*\
73                            Class directionInfo Declaration
74 \*---------------------------------------------------------------------------*/
76 class directionInfo
78     // Private data
80         // Either mesh edge or face point
81         label index_;
83         // Local n axis
84         vector n_;
87     // Private Member Functions
89         //- Find edge among edgeLabels that uses v0 and v1
90         static label findEdge
91         (
92             const primitiveMesh& mesh,
93             const labelList& edgeLabels,
94             const label v1,
95             const label v0
96         );
98         //- Return 'lowest' of a,b in face of size.
99         static label lowest
100         (
101             const label size,
102             const label a,
103             const label b
104         );
106 public:
108     // Static Functions
110         //- Given edge on hex cell find corresponding edge on face. Is either
111         //  index in face or -1 (cut tangential to face). Public since is
112         //  needed to fill in seed faces in meshWave.
113         static label edgeToFaceIndex
114         (
115             const primitiveMesh& mesh,
116             const label cellI,
117             const label faceI,
118             const label edgeI
119         );
121     // Constructors
123         //- Construct null
124         inline directionInfo();
126         //- Construct from components
127         inline directionInfo
128         (
129             const label,
130             const vector& n
131         );
133         //- Construct as copy
134         inline directionInfo(const directionInfo&);
137     // Member Functions
139         // Access
141             inline label index() const
142             {
143                 return index_;
144             }
146             inline const vector& n() const
147             {
148                 return n_;
149             }
151         // Needed by FaceCellWave
153             //- Check whether origin has been changed at all or
154             //  still contains original (invalid) value.
155             template<class TrackingData>
156             inline bool valid(TrackingData& td) const;
158             //- Check for identical geometrical data. Used for cyclics checking.
159             template<class TrackingData>
160             inline bool sameGeometry
161             (
162                 const polyMesh&,
163                 const directionInfo&,
164                 const scalar,
165                 TrackingData& td
166             ) const;
168             //- Convert any absolute coordinates into relative to (patch)face
169             //  centre
170             template<class TrackingData>
171             inline void leaveDomain
172             (
173                 const polyMesh&,
174                 const polyPatch&,
175                 const label patchFaceI,
176                 const point& faceCentre,
177                 TrackingData& td
178             );
180             //- Reverse of leaveDomain
181             template<class TrackingData>
182             inline void enterDomain
183             (
184                 const polyMesh&,
185                 const polyPatch&,
186                 const label patchFaceI,
187                 const point& faceCentre,
188                 TrackingData& td
189             );
191             //- Apply rotation matrix to any coordinates
192             template<class TrackingData>
193             inline void transform
194             (
195                 const polyMesh&,
196                 const tensor&,
197                 TrackingData& td
198             );
200             //- Influence of neighbouring face.
201             template<class TrackingData>
202             inline bool updateCell
203             (
204                 const polyMesh&,
205                 const label thisCellI,
206                 const label neighbourFaceI,
207                 const directionInfo& neighbourInfo,
208                 const scalar tol,
209                 TrackingData& td
210             );
212             //- Influence of neighbouring cell.
213             template<class TrackingData>
214             inline bool updateFace
215             (
216                 const polyMesh&,
217                 const label thisFaceI,
218                 const label neighbourCellI,
219                 const directionInfo& neighbourInfo,
220                 const scalar tol,
221                 TrackingData& td
222             );
224             //- Influence of different value on same face.
225             template<class TrackingData>
226             inline bool updateFace
227             (
228                 const polyMesh&,
229                 const label thisFaceI,
230                 const directionInfo& neighbourInfo,
231                 const scalar tol,
232                 TrackingData& td
233             );
235             //- Same (like operator==)
236             template<class TrackingData>
237             inline bool equal(const directionInfo&, TrackingData& td) const;
239     // Member Operators
241         // Needed for List IO
242         inline bool operator==(const directionInfo&) const;
244         inline bool operator!=(const directionInfo&) const;
247     // IOstream Operators
249         friend Ostream& operator<<(Ostream&, const directionInfo&);
250         friend Istream& operator>>(Istream&, directionInfo&);
254 //- Data associated with directionInfo type are contiguous
255 template<>
256 inline bool contiguous<directionInfo>()
258     return true;
262 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 } // End namespace Foam
266 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268 #include "directionInfoI.H"
270 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 #endif
274 // ************************************************************************* //