directionInfo.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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Class
25  Foam::directionInfo
26 
27 Description
28  Holds direction in which to split cell (in fact a local coordinate axes).
29  Information is a label and a direction.
30 
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
38 
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)
41 
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.
47 
48 SourceFiles
49  directionInfoI.H
50  directionInfo.C
51 
52 \*---------------------------------------------------------------------------*/
53 
54 #ifndef directionInfo_H
55 #define directionInfo_H
56 
57 #include "point.H"
58 #include "labelList.H"
59 #include "tensor.H"
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 namespace Foam
64 {
65 class polyPatch;
66 class polyMesh;
67 class primitiveMesh;
68 class edge;
69 class face;
70 class polyMesh;
71 
72 
73 // Forward declaration of friend functions and operators
74 
75 class directionInfo;
76 
77 Istream& operator>>(Istream&, directionInfo&);
78 Ostream& operator<<(Ostream&, const directionInfo&);
79 
80 
81 /*---------------------------------------------------------------------------*\
82  Class directionInfo Declaration
83 \*---------------------------------------------------------------------------*/
84 
85 class directionInfo
86 {
87  // Private data
88 
89  // Either mesh edge or face point
90  label index_;
91 
92  // Local n axis
93  vector n_;
94 
95 
96  // Private Member Functions
97 
98  //- Find edge among edgeLabels that uses v0 and v1
99  static label findEdge
100  (
101  const primitiveMesh& mesh,
102  const labelList& edgeLabels,
103  const label v1,
104  const label v0
105  );
106 
107  //- Return 'lowest' of a,b in face of size.
108  static label lowest
109  (
110  const label size,
111  const label a,
112  const label b
113  );
114 
115 public:
116 
117  // Static Functions
118 
119  //- Given edge on hex cell find corresponding edge on face. Is either
120  // index in face or -1 (cut tangential to face). Public since is
121  // needed to fill in seed faces in meshWave.
122  static label edgeToFaceIndex
123  (
124  const primitiveMesh& mesh,
125  const label celli,
126  const label facei,
127  const label edgeI
128  );
129 
130  // Constructors
131 
132  //- Construct null
133  inline directionInfo();
134 
135  //- Construct from components
136  inline directionInfo
137  (
138  const label,
139  const vector& n
140  );
141 
142  //- Construct as copy
143  inline directionInfo(const directionInfo&);
144 
145 
146  // Member Functions
147 
148  // Access
150  inline label index() const
151  {
152  return index_;
153  }
155  inline const vector& n() const
156  {
157  return n_;
158  }
159 
160  // Needed by FaceCellWave
161 
162  //- Check whether origin has been changed at all or
163  // still contains original (invalid) value.
164  template<class TrackingData>
165  inline bool valid(TrackingData& td) const;
166 
167  //- Check for identical geometrical data. Used for cyclics checking.
168  template<class TrackingData>
169  inline bool sameGeometry
170  (
171  const polyMesh&,
172  const directionInfo&,
173  const scalar,
174  TrackingData& td
175  ) const;
176 
177  //- Convert any absolute coordinates into relative to (patch)face
178  // centre
179  template<class TrackingData>
180  inline void leaveDomain
181  (
182  const polyMesh&,
183  const polyPatch&,
184  const label patchFacei,
185  const point& faceCentre,
186  TrackingData& td
187  );
188 
189  //- Reverse of leaveDomain
190  template<class TrackingData>
191  inline void enterDomain
192  (
193  const polyMesh&,
194  const polyPatch&,
195  const label patchFacei,
196  const point& faceCentre,
197  TrackingData& td
198  );
199 
200  //- Apply rotation matrix to any coordinates
201  template<class TrackingData>
202  inline void transform
203  (
204  const polyMesh&,
205  const tensor&,
206  TrackingData& td
207  );
208 
209  //- Influence of neighbouring face.
210  template<class TrackingData>
211  inline bool updateCell
212  (
213  const polyMesh&,
214  const label thisCelli,
215  const label neighbourFacei,
216  const directionInfo& neighbourInfo,
217  const scalar tol,
218  TrackingData& td
219  );
220 
221  //- Influence of neighbouring cell.
222  template<class TrackingData>
223  inline bool updateFace
224  (
225  const polyMesh&,
226  const label thisFacei,
227  const label neighbourCelli,
228  const directionInfo& neighbourInfo,
229  const scalar tol,
230  TrackingData& td
231  );
232 
233  //- Influence of different value on same face.
234  template<class TrackingData>
235  inline bool updateFace
236  (
237  const polyMesh&,
238  const label thisFacei,
239  const directionInfo& neighbourInfo,
240  const scalar tol,
241  TrackingData& td
242  );
243 
244  //- Same (like operator==)
245  template<class TrackingData>
246  inline bool equal(const directionInfo&, TrackingData& td) const;
247 
248  // Member Operators
249 
250  // Needed for List IO
251  inline bool operator==(const directionInfo&) const;
252 
253  inline bool operator!=(const directionInfo&) const;
254 
255 
256  // IOstream Operators
257 
258  friend Ostream& operator<<(Ostream&, const directionInfo&);
260 };
261 
262 
263 //- Data associated with directionInfo type are contiguous
264 template<>
265 inline bool contiguous<directionInfo>()
266 {
267  return true;
268 }
269 
270 
271 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 
273 } // End namespace Foam
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 #include "directionInfoI.H"
278 
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280 
281 #endif
282 
283 // ************************************************************************* //
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:59
bool equal(const directionInfo &, TrackingData &td) const
Same (like operator==)
const vector & n() const
label index() const
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
friend Ostream & operator<<(Ostream &, const directionInfo &)
void leaveDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Convert any absolute coordinates into relative to (patch)face.
bool updateCell(const polyMesh &, const label thisCelli, const label neighbourFacei, const directionInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.
bool contiguous< directionInfo >()
Data associated with directionInfo type are contiguous.
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
dynamicFvMesh & mesh
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
void enterDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Reverse of leaveDomain.
bool operator!=(const directionInfo &) const
Istream & operator>>(Istream &, directionInfo &)
bool operator==(const directionInfo &) const
static label edgeToFaceIndex(const primitiveMesh &mesh, const label celli, const label facei, const label edgeI)
Given edge on hex cell find corresponding edge on face. Is either.
Definition: directionInfo.C:93
friend Istream & operator>>(Istream &, directionInfo &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Holds direction in which to split cell (in fact a local coordinate axes). Information is a label and ...
Definition: directionInfo.H:84
directionInfo()
Construct null.
void transform(const polyMesh &, const tensor &, TrackingData &td)
Apply rotation matrix to any coordinates.
bool updateFace(const polyMesh &, const label thisFacei, const label neighbourCelli, const directionInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring cell.
Ostream & operator<<(Ostream &, const ensightPart &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
bool sameGeometry(const polyMesh &, const directionInfo &, const scalar, TrackingData &td) const
Check for identical geometrical data. Used for cyclics checking.
Namespace for OpenFOAM.