patchFaceOrientationI.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) 2013-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 \*---------------------------------------------------------------------------*/
25 
26 #include "polyMesh.H"
27 #include "transform.H"
28 #include "orientedSurface.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
33 :
34  flipStatus_(orientedSurface::UNVISITED)
35 {}
36 
37 
39 (
40  const label flipStatus
41 )
42 :
43  flipStatus_(flipStatus)
44 {}
45 
46 
47 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
48 
50 {
51  return flipStatus_;
52 }
53 
54 
56 {
57  if (flipStatus_ == orientedSurface::NOFLIP)
58  {
59  flipStatus_ = orientedSurface::FLIP;
60  }
61  else if (flipStatus_ == orientedSurface::FLIP)
62  {
63  flipStatus_ = orientedSurface::NOFLIP;
64  }
65 }
66 
67 
68 template<class TrackingData>
69 inline bool Foam::patchFaceOrientation::valid(TrackingData& td) const
70 {
71  return flipStatus_ != orientedSurface::UNVISITED;
72 }
73 
74 
75 template<class TrackingData>
77 (
78  const polyMesh& mesh,
79  const indirectPrimitivePatch& patch,
80  const tensor& rotTensor,
81  const scalar tol,
82  TrackingData& td
83 )
84 {}
85 
86 
87 template<class TrackingData>
89 (
90  const polyMesh& mesh,
91  const indirectPrimitivePatch& patch,
92  const label edgeI,
93  const label facei,
94  const patchFaceOrientation& faceInfo,
95  const scalar tol,
96  TrackingData& td
97 )
98 {
99  if (valid(td))
100  {
101  return false;
102  }
103 
104  const face& f = patch.localFaces()[facei];
105  const edge& e = patch.edges()[edgeI];
106 
107  //Pout<< "Updating edge:" << edgeI << " verts:" << e << nl
108  // << " start:" << patch.localPoints()[e[0]] << nl
109  // << " end:" << patch.localPoints()[e[1]] << endl;
110 
111  patchFaceOrientation consistentInfo(faceInfo);
112 
113  // Check how edge relates to face
114  if (f.edgeDirection(e) < 0)
115  {
116  // Create flipped version of faceInfo
117  consistentInfo.flip();
118  }
119 
120  operator=(consistentInfo);
121  return true;
122 }
123 
124 
125 template<class TrackingData>
127 (
128  const polyMesh& mesh,
129  const indirectPrimitivePatch& patch,
130  const patchFaceOrientation& edgeInfo,
131  const bool sameOrientation,
132  const scalar tol,
133  TrackingData& td
134 )
135 {
136  if (valid(td))
137  {
138  return false;
139  }
140 
141  // Create (flipped/unflipped) version of edgeInfo
142  patchFaceOrientation consistentInfo(edgeInfo);
143 
144  if (!sameOrientation)
145  {
146  consistentInfo.flip();
147  }
148 
149  operator=(consistentInfo);
150  return true;
151 }
152 
153 
154 template<class TrackingData>
156 (
157  const polyMesh& mesh,
158  const indirectPrimitivePatch& patch,
159  const label facei,
160  const label edgeI,
161  const patchFaceOrientation& edgeInfo,
162  const scalar tol,
163  TrackingData& td
164 )
165 {
166  if (valid(td))
167  {
168  return false;
169  }
170 
171  // Transfer flip to face
172  const face& f = patch.localFaces()[facei];
173  const edge& e = patch.edges()[edgeI];
174 
175 
176  //Pout<< "Updating face:" << facei << nl
177  // << " verts:" << f << nl
178  // << " with edge:" << edgeI << nl
179  // << " start:" << patch.localPoints()[e[0]] << nl
180  // << " end:" << patch.localPoints()[e[1]] << endl;
181 
182 
183  // Create (flipped/unflipped) version of edgeInfo
184  patchFaceOrientation consistentInfo(edgeInfo);
185 
186  if (f.edgeDirection(e) > 0)
187  {
188  consistentInfo.flip();
189  }
190 
191  operator=(consistentInfo);
192  return true;
193 }
194 
195 
196 template<class TrackingData>
198 (
199  const patchFaceOrientation& rhs,
200  TrackingData& td
201 ) const
202 {
203  return operator==(rhs);
204 }
205 
206 
207 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
208 
209 inline bool Foam::patchFaceOrientation::operator==
210 (
211  const Foam::patchFaceOrientation& rhs
212 ) const
213 {
214  return flipStatus() == rhs.flipStatus();
215 }
216 
217 
218 inline bool Foam::patchFaceOrientation::operator!=
219 (
220  const Foam::patchFaceOrientation& rhs
221 ) const
222 {
223  return !(*this == rhs);
224 }
225 
226 
227 // ************************************************************************* //
patchFaceOrientation()
Construct null.
int edgeDirection(const edge &) const
Return the edge direction on the face.
Definition: face.C:779
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 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.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const double e
Elementary charge.
Definition: doubleFloat.H:78
bool operator==(const patchFaceOrientation &) const
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Given point flip all faces such that normals point in same direction.
Transport of orientation for use in PatchEdgeFaceWave.
bool equal(const patchFaceOrientation &, TrackingData &) const
Same (like operator==)
A list of faces which address into the list of points.
label flipStatus() const
Orientation.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
3D tensor transformation operations.
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.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
void flip()
Reverse orientation.
labelList f(nPoints)
void transform(const polyMesh &mesh, const indirectPrimitivePatch &patch, const tensor &rotTensor, const scalar tol, TrackingData &td)
Apply rotation matrix.
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74