1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
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
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/>.
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
52 \*---------------------------------------------------------------------------*/
54 #ifndef directionInfo_H
55 #define directionInfo_H
58 #include "labelList.H"
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 /*---------------------------------------------------------------------------*\
73 Class directionInfo Declaration
74 \*---------------------------------------------------------------------------*/
80 // Either mesh edge or face point
87 // Private Member Functions
89 //- Find edge among edgeLabels that uses v0 and v1
92 const primitiveMesh& mesh,
93 const labelList& edgeLabels,
98 //- Return 'lowest' of a,b in face of size.
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
115 const primitiveMesh& mesh,
124 inline directionInfo();
126 //- Construct from components
133 //- Construct as copy
134 inline directionInfo(const directionInfo&);
141 inline label index() const
146 inline const vector& n() const
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
163 const directionInfo&,
168 //- Convert any absolute coordinates into relative to (patch)face
170 template<class TrackingData>
171 inline void leaveDomain
175 const label patchFaceI,
176 const point& faceCentre,
180 //- Reverse of leaveDomain
181 template<class TrackingData>
182 inline void enterDomain
186 const label patchFaceI,
187 const point& faceCentre,
191 //- Apply rotation matrix to any coordinates
192 template<class TrackingData>
193 inline void transform
200 //- Influence of neighbouring face.
201 template<class TrackingData>
202 inline bool updateCell
205 const label thisCellI,
206 const label neighbourFaceI,
207 const directionInfo& neighbourInfo,
212 //- Influence of neighbouring cell.
213 template<class TrackingData>
214 inline bool updateFace
217 const label thisFaceI,
218 const label neighbourCellI,
219 const directionInfo& neighbourInfo,
224 //- Influence of different value on same face.
225 template<class TrackingData>
226 inline bool updateFace
229 const label thisFaceI,
230 const directionInfo& neighbourInfo,
235 //- Same (like operator==)
236 template<class TrackingData>
237 inline bool equal(const directionInfo&, TrackingData& td) const;
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
256 inline bool contiguous<directionInfo>()
262 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 } // End namespace Foam
266 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268 #include "directionInfoI.H"
270 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274 // ************************************************************************* //