patchFaceOrientationI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2021 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 \*---------------------------------------------------------------------------*/
25 
26 #include "orientedSurface.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
31 :
32  flipStatus_(orientedSurface::UNVISITED)
33 {}
34 
35 
37 (
38  const label flipStatus
39 )
40 :
41  flipStatus_(flipStatus)
42 {}
43 
44 
45 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
46 
48 {
49  return flipStatus_;
50 }
51 
52 
54 {
55  if (flipStatus_ == orientedSurface::NOFLIP)
56  {
57  flipStatus_ = orientedSurface::FLIP;
58  }
59  else if (flipStatus_ == orientedSurface::FLIP)
60  {
61  flipStatus_ = orientedSurface::NOFLIP;
62  }
63 }
64 
65 
66 template<class TrackingData>
67 inline bool Foam::patchFaceOrientation::valid(TrackingData& td) const
68 {
69  return flipStatus_ != orientedSurface::UNVISITED;
70 }
71 
72 
73 template<class TrackingData>
75 (
76  const polyMesh& mesh,
77  const indirectPrimitivePatch& patch,
78  const tensor& rotTensor,
79  const scalar tol,
80  TrackingData& td
81 )
82 {}
83 
84 
85 template<class TrackingData>
87 (
88  const polyMesh& mesh,
89  const indirectPrimitivePatch& patch,
90  const label edgei,
91  const label facei,
92  const patchFaceOrientation& faceInfo,
93  const scalar tol,
94  TrackingData& td
95 )
96 {
97  if (valid(td))
98  {
99  return false;
100  }
101 
102  const face& f = patch.localFaces()[facei];
103  const edge& e = patch.edges()[edgei];
104 
105  patchFaceOrientation consistentInfo(faceInfo);
106 
107  // Check how edge relates to face
108  if (f.edgeDirection(e) < 0)
109  {
110  // Create flipped version of faceInfo
111  consistentInfo.flip();
112  }
113 
114  operator=(consistentInfo);
115 
116  return true;
117 }
118 
119 
120 template<class TrackingData>
122 (
123  const polyMesh& mesh,
124  const indirectPrimitivePatch& patch,
125  const patchFaceOrientation& edgeinfo,
126  const bool sameOrientation,
127  const scalar tol,
128  TrackingData& td
129 )
130 {
131  if (valid(td))
132  {
133  return false;
134  }
135 
136  // Create (flipped/unflipped) version of edgeinfo
137  patchFaceOrientation consistentInfo(edgeinfo);
138 
139  if (!sameOrientation)
140  {
141  consistentInfo.flip();
142  }
143 
144  operator=(consistentInfo);
145 
146  return true;
147 }
148 
149 
150 template<class TrackingData>
152 (
153  const polyMesh& mesh,
154  const indirectPrimitivePatch& patch,
155  const label facei,
156  const label edgei,
157  const patchFaceOrientation& edgeinfo,
158  const scalar tol,
159  TrackingData& td
160 )
161 {
162  if (valid(td))
163  {
164  return false;
165  }
166 
167  // Transfer flip to face
168  const face& f = patch.localFaces()[facei];
169  const edge& e = patch.edges()[edgei];
170 
171  // Create (flipped/unflipped) version of edgeinfo
172  patchFaceOrientation consistentInfo(edgeinfo);
173 
174  if (f.edgeDirection(e) > 0)
175  {
176  consistentInfo.flip();
177  }
178 
179  operator=(consistentInfo);
180 
181  return true;
182 }
183 
184 
185 template<class TrackingData>
187 (
188  const patchFaceOrientation& rhs,
189  TrackingData& td
190 ) const
191 {
192  return operator==(rhs);
193 }
194 
195 
196 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
197 
198 inline bool Foam::patchFaceOrientation::operator==
199 (
200  const Foam::patchFaceOrientation& rhs
201 ) const
202 {
203  return flipStatus() == rhs.flipStatus();
204 }
205 
206 
207 inline bool Foam::patchFaceOrientation::operator!=
208 (
209  const Foam::patchFaceOrientation& rhs
210 ) const
211 {
212  return !(*this == rhs);
213 }
214 
215 
216 // ************************************************************************* //
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Given point flip all faces such that normals point in same direction.
Transport of orientation for use in PatchEdgeFaceWave.
void flip()
Reverse orientation.
bool updateFace(const polyMesh &mesh, const indirectPrimitivePatch &patch, const label facei, const label edgei, const patchFaceOrientation &edgeinfo, const scalar tol, TrackingData &td)
Influence of edge on face.
bool updateEdge(const polyMesh &mesh, const indirectPrimitivePatch &patch, const label edgei, const label facei, const patchFaceOrientation &faceInfo, const scalar tol, TrackingData &td)
Influence of face on edge.
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
label flipStatus() const
Orientation.
bool equal(const patchFaceOrientation &, TrackingData &) const
Same (like operator==)
void transform(const polyMesh &mesh, const indirectPrimitivePatch &patch, const tensor &rotTensor, const scalar tol, TrackingData &td)
Apply rotation matrix.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
bool valid(const PtrList< ModelType > &l)
const doubleScalar e
Definition: doubleScalar.H:106
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
labelList f(nPoints)