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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::directionInfo
29 
30 Description
31  Holds direction in which to split cell (in fact a local coordinate axes).
32  Information is a label and a direction.
33 
34  The direction is the normal
35  direction to cut in. The label's meaning depends on whether the info
36  is on a cell or on a face:
37  - in cell: edge that is being cut. (determines for hex how cut is)
38  - in face: local face point that is being cut or -1.
39  -# (-1) : cut is tangential to plane
40  -# (>= 0): edge fp..fp+1 is cut
41 
42  (has to be facepoint, not vertex since vertex not valid across
43  processors whereas f[0] should correspond to f[0] on other side)
44 
45  The rule is that if the label is set (-1 or higher) it is used
46  (topological information only), otherwise the vector is used. This makes
47  sure that we use topological information as much as possible and so a
48  hex mesh is cut purely topologically. All other shapes are cut
49  geometrically.
50 
51 SourceFiles
52  directionInfoI.H
53  directionInfo.C
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #ifndef directionInfo_H
58 #define directionInfo_H
59 
60 #include "point.H"
61 #include "labelList.H"
62 #include "tensor.H"
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 namespace Foam
67 {
68 
69 // Forward Declarations
70 class polyPatch;
71 class polyMesh;
72 class primitiveMesh;
73 class edge;
74 class face;
75 class polyMesh;
76 class directionInfo;
77 
78 Istream& operator>>(Istream&, directionInfo&);
79 Ostream& operator<<(Ostream&, const directionInfo&);
80 
81 
82 /*---------------------------------------------------------------------------*\
83  Class directionInfo Declaration
84 \*---------------------------------------------------------------------------*/
85 
86 class directionInfo
87 {
88  // Private Data
89 
90  // Either mesh edge or face point
91  label index_;
92 
93  // Local n axis
94  vector n_;
95 
96 
97  // Private Member Functions
98 
99  //- Find edge among edgeLabels that uses v0 and v1
100  static label findEdge
101  (
102  const primitiveMesh& mesh,
103  const labelList& edgeLabels,
104  const label v1,
105  const label v0
106  );
107 
108  //- Return 'lowest' of a,b in face of size.
109  static label lowest
110  (
111  const label size,
112  const label a,
113  const label b
114  );
115 
116 public:
117 
118  // Static Functions
119 
120  //- Given edge on hex cell find corresponding edge on face. Is either
121  // index in face or -1 (cut tangential to face). Public since is
122  // needed to fill in seed faces in meshWave.
123  static label edgeToFaceIndex
124  (
125  const primitiveMesh& mesh,
126  const label celli,
127  const label facei,
128  const label edgeI
129  );
130 
131  // Constructors
132 
133  //- Construct null
134  inline directionInfo();
135 
136  //- Construct from components
137  inline directionInfo
138  (
139  const label,
140  const vector& n
141  );
142 
143  //- Construct as copy
144  inline directionInfo(const directionInfo&);
145 
146 
147  // Member Functions
148 
149  // Access
150 
151  inline label index() const
152  {
153  return index_;
154  }
155 
156  inline const vector& n() const
157  {
158  return n_;
159  }
160 
161  // Needed by FaceCellWave
162 
163  //- Check whether origin has been changed at all or
164  // still contains original (invalid) value.
165  template<class TrackingData>
166  inline bool valid(TrackingData& td) const;
167 
168  //- Check for identical geometrical data. Used for cyclics checking.
169  template<class TrackingData>
170  inline bool sameGeometry
171  (
172  const polyMesh&,
173  const directionInfo&,
174  const scalar,
175  TrackingData& td
176  ) const;
177 
178  //- Convert any absolute coordinates into relative to (patch)face
179  // centre
180  template<class TrackingData>
181  inline void leaveDomain
182  (
183  const polyMesh&,
184  const polyPatch&,
185  const label patchFacei,
186  const point& faceCentre,
187  TrackingData& td
188  );
189 
190  //- Reverse of leaveDomain
191  template<class TrackingData>
192  inline void enterDomain
193  (
194  const polyMesh&,
195  const polyPatch&,
196  const label patchFacei,
197  const point& faceCentre,
198  TrackingData& td
199  );
200 
201  //- Apply rotation matrix to any coordinates
202  template<class TrackingData>
203  inline void transform
204  (
205  const polyMesh&,
206  const tensor&,
207  TrackingData& td
208  );
209 
210  //- Influence of neighbouring face.
211  template<class TrackingData>
212  inline bool updateCell
213  (
214  const polyMesh&,
215  const label thisCelli,
216  const label neighbourFacei,
217  const directionInfo& neighbourInfo,
218  const scalar tol,
219  TrackingData& td
220  );
221 
222  //- Influence of neighbouring cell.
223  template<class TrackingData>
224  inline bool updateFace
225  (
226  const polyMesh&,
227  const label thisFacei,
228  const label neighbourCelli,
229  const directionInfo& neighbourInfo,
230  const scalar tol,
231  TrackingData& td
232  );
233 
234  //- Influence of different value on same face.
235  template<class TrackingData>
236  inline bool updateFace
237  (
238  const polyMesh&,
239  const label thisFacei,
240  const directionInfo& neighbourInfo,
241  const scalar tol,
242  TrackingData& td
243  );
244 
245  //- Same (like operator==)
246  template<class TrackingData>
247  inline bool equal(const directionInfo&, TrackingData& td) const;
248 
249  // Member Operators
250 
251  // Needed for List IO
252  inline bool operator==(const directionInfo&) const;
253 
254  inline bool operator!=(const directionInfo&) const;
255 
256 
257  // IOstream Operators
258 
259  friend Ostream& operator<<(Ostream&, const directionInfo&);
261 };
262 
263 
264 // * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
265 
266 //- Contiguous data for directionInfo
267 template<> struct is_contiguous<directionInfo> : std::true_type {};
268 
269 
270 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
271 
272 } // End namespace Foam
273 
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 
276 #include "directionInfoI.H"
277 
278 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279 
280 #endif
281 
282 // ************************************************************************* //
Foam::directionInfo::edgeToFaceIndex
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:95
Foam::Tensor< scalar >
Foam::directionInfo::index
label index() const
Definition: directionInfo.H:150
Foam::directionInfo::enterDomain
void enterDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Reverse of leaveDomain.
Definition: directionInfoI.H:105
point.H
Foam::directionInfo::equal
bool equal(const directionInfo &, TrackingData &td) const
Same (like operator==)
Definition: directionInfoI.H:292
Foam::directionInfo::operator!=
bool operator!=(const directionInfo &) const
Definition: directionInfoI.H:310
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:217
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:228
directionInfoI.H
tensor.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::directionInfo::valid
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
Definition: directionInfoI.H:65
Foam::directionInfo::transform
void transform(const polyMesh &, const tensor &, TrackingData &td)
Apply rotation matrix to any coordinates.
Definition: directionInfoI.H:125
Foam::label
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:62
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
labelList.H
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::directionInfo
Holds direction in which to split cell (in fact a local coordinate axes). Information is a label and ...
Definition: directionInfo.H:85
Foam::directionInfo::operator<<
friend Ostream & operator<<(Ostream &, const directionInfo &)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::directionInfo::operator>>
friend Istream & operator>>(Istream &, directionInfo &)
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::directionInfo::n
const vector & n() const
Definition: directionInfo.H:155
Foam::directionInfo::operator==
bool operator==(const directionInfo &) const
Definition: directionInfoI.H:303
Foam::Vector< scalar >
Foam::List< label >
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:89
Foam::directionInfo::directionInfo
directionInfo()
Construct null.
Definition: directionInfoI.H:35
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &)
Definition: boundaryPatch.C:102
Foam::directionInfo::sameGeometry
bool sameGeometry(const polyMesh &, const directionInfo &, const scalar, TrackingData &td) const
Check for identical geometrical data. Used for cyclics checking.
Definition: directionInfoI.H:74
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:136
Foam::is_contiguous
A template class to specify that a data type can be considered as being contiguous in memory.
Definition: contiguous.H:75
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:78