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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
28
29#include "polyMesh.H"
30#include "transform.H"
31#include "orientedSurface.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
36:
37 flipStatus_(orientedSurface::UNVISITED)
38{}
39
40
42(
43 const label flipStatus
44)
45:
46 flipStatus_(flipStatus)
47{}
48
49
50// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
51
53{
54 if (flipStatus_ == orientedSurface::NOFLIP)
55 {
56 flipStatus_ = orientedSurface::FLIP;
57 }
58 else if (flipStatus_ == orientedSurface::FLIP)
59 {
60 flipStatus_ = orientedSurface::NOFLIP;
61 }
62}
63
64
65template<class TrackingData>
66inline bool Foam::patchFaceOrientation::valid(TrackingData& td) const
67{
68 return flipStatus_ != orientedSurface::UNVISITED;
69}
70
71
72template<class TrackingData>
74(
75 const polyMesh& mesh,
76 const indirectPrimitivePatch& patch,
77 const tensor& rotTensor,
78 const scalar tol,
79 TrackingData& td
80)
81{}
82
83
84template<class TrackingData>
86(
87 const polyMesh& mesh,
88 const indirectPrimitivePatch& patch,
89 const label edgeI,
90 const label facei,
91 const patchFaceOrientation& faceInfo,
92 const scalar tol,
93 TrackingData& td
94)
95{
96 if (valid(td))
97 {
98 return false;
99 }
100
101 const face& f = patch.localFaces()[facei];
102 const edge& e = patch.edges()[edgeI];
103
104 //Pout<< "Updating edge:" << edgeI << " verts:" << e << nl
105 // << " start:" << patch.localPoints()[e[0]] << nl
106 // << " end:" << patch.localPoints()[e[1]] << endl;
107
108 patchFaceOrientation consistentInfo(faceInfo);
109
110 // Check how edge relates to face
111 if (f.edgeDirection(e) < 0)
112 {
113 // Create flipped version of faceInfo
114 consistentInfo.flip();
115 }
116
117 operator=(consistentInfo);
118 return true;
119}
120
121
122template<class TrackingData>
124(
125 const polyMesh& mesh,
126 const indirectPrimitivePatch& patch,
127 const patchFaceOrientation& edgeInfo,
128 const bool sameOrientation,
129 const scalar tol,
130 TrackingData& td
131)
132{
133 if (valid(td))
134 {
135 return false;
136 }
137
138 // Create (flipped/unflipped) version of edgeInfo
139 patchFaceOrientation consistentInfo(edgeInfo);
140
141 if (!sameOrientation)
142 {
143 consistentInfo.flip();
144 }
145
146 operator=(consistentInfo);
147 return true;
148}
149
150
151template<class TrackingData>
153(
154 const polyMesh& mesh,
155 const indirectPrimitivePatch& patch,
156 const label facei,
157 const label edgeI,
158 const patchFaceOrientation& edgeInfo,
159 const scalar tol,
160 TrackingData& td
161)
162{
163 if (valid(td))
164 {
165 return false;
166 }
167
168 // Transfer flip to face
169 const face& f = patch.localFaces()[facei];
170 const edge& e = patch.edges()[edgeI];
171
172
173 //Pout<< "Updating face:" << facei << nl
174 // << " verts:" << f << nl
175 // << " with edge:" << edgeI << nl
176 // << " start:" << patch.localPoints()[e[0]] << nl
177 // << " end:" << patch.localPoints()[e[1]] << endl;
178
179
180 // Create (flipped/unflipped) version of edgeInfo
181 patchFaceOrientation consistentInfo(edgeInfo);
182
183 if (f.edgeDirection(e) > 0)
184 {
185 consistentInfo.flip();
186 }
187
188 operator=(consistentInfo);
189 return true;
190}
191
192
193template<class TrackingData>
195(
196 const patchFaceOrientation& rhs,
197 TrackingData& td
198) const
199{
200 return operator==(rhs);
201}
202
203
204// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
205
207(
208 const patchFaceOrientation& rhs
209) const
210{
211 return flipStatus_ == rhs.flipStatus_;
212}
213
214
216(
217 const patchFaceOrientation& rhs
218) const
219{
220 return !(*this == rhs);
221}
222
223
224// ************************************************************************* //
A list of faces which address into the list of points.
bool valid() const
True if all internal ids are non-negative.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Default transformation behaviour.
Given point flip all faces such that normals point in same direction.
Transport of orientation for use in PatchEdgeFaceWave.
void flip()
Reverse the 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.
patchFaceOrientation()
Default construct.
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 equal(const patchFaceOrientation &, TrackingData &) const
Test for equality, with TrackingData.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
dynamicFvMesh & mesh
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
labelList f(nPoints)
volScalarField & e
Definition: createFields.H:11
3D tensor transformation operations.