Loading...
Searching...
No Matches
directionInfo.C
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) 2020-2022 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 "directionInfo.H"
30#include "polyMesh.H"
31
32// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33
34// Find edge among edgeLabels that uses v0 and v1
35Foam::label Foam::directionInfo::findEdge
36(
37 const primitiveMesh& mesh,
38 const labelList& edgeLabels,
39 const label v1,
40 const label v0
41)
42{
43 forAll(edgeLabels, edgeLabelI)
44 {
45 label edgeI = edgeLabels[edgeLabelI];
46
47 if (mesh.edges()[edgeI] == edge(v0, v1))
48 {
49 return edgeI;
50 }
51 }
52
54 << "Cannot find an edge among " << edgeLabels << endl
55 << "that uses vertices " << v0
56 << " and " << v1
57 << abort(FatalError);
58
59 return -1;
60}
61
62
63Foam::label Foam::directionInfo::lowest
64(
65 const label size,
66 const label a,
67 const label b
68)
69{
70 // Get next point
71 label a1 = (a + 1) % size;
72
73 if (a1 == b)
74 {
75 return a;
76 }
77 else
78 {
79 label b1 = (b + 1) % size;
80
81 if (b1 != a)
82 {
84 << "Problem : a:" << a << " b:" << b << " size:" << size
85 << abort(FatalError);
86 }
87
88 return b;
89 }
90}
91
92
93// Have edge on hex cell. Find corresponding edge on face. Return -1 if none
94// found.
96(
97 const primitiveMesh& mesh,
98 const label celli,
99 const label facei,
100 const label edgeI
101)
102{
103 if ((edgeI < 0) || (edgeI >= mesh.nEdges()))
104 {
106 << "Illegal edge label:" << edgeI
107 << " when projecting cut edge from cell " << celli
108 << " to face " << facei
109 << abort(FatalError);
110 }
111
112 const edge& e = mesh.edges()[edgeI];
113
114 const face& f = mesh.faces()[facei];
115
116 // edgeI is either
117 // - in facei. Convert into index in face.
118 // - connected (but not in) to face. Return -1.
119 // - in face opposite facei. Convert into index in face.
120
121 label fpA = f.find(e.start());
122 label fpB = f.find(e.end());
123
124 if (fpA != -1)
125 {
126 if (fpB != -1)
127 {
128 return lowest(f.size(), fpA, fpB);
129 }
130 else
131 {
132 // e.start() in face, e.end() not
133 return -1;
134 }
135 }
136 else
137 {
138 if (fpB != -1)
139 {
140 // e.end() in face, e.start() not
141 return -1;
142 }
143 else
144 {
145 // Both not in face.
146 // e is on opposite face. Determine corresponding edge on this face:
147 // - determine two faces using edge (one is the opposite face,
148 // one is 'side' face
149 // - walk on both these faces to opposite edge
150 // - check if this opposite edge is on facei
151
152 label f0I, f1I;
153
154 meshTools::getEdgeFaces(mesh, celli, edgeI, f0I, f1I);
155
156 // Walk to opposite edge on face f0
157 label edge0I =
158 meshTools::walkFace(mesh, f0I, edgeI, e.start(), 2);
159
160 // Check if edge on facei.
161
162 const edge& e0 = mesh.edges()[edge0I];
163
164 fpA = f.find(e0.start());
165 fpB = f.find(e0.end());
166
167 if ((fpA != -1) && (fpB != -1))
168 {
169 return lowest(f.size(), fpA, fpB);
170 }
171
172 // Face0 is doesn't have an edge on facei (so must be the opposite
173 // face) so try face1.
174
175 // Walk to opposite edge on face f1
176 label edge1I =
177 meshTools::walkFace(mesh, f1I, edgeI, e.start(), 2);
178
179 // Check if edge on facei.
180 const edge& e1 = mesh.edges()[edge1I];
181
182 fpA = f.find(e1.start());
183 fpB = f.find(e1.end());
184
185 if ((fpA != -1) && (fpB != -1))
186 {
187 return lowest(f.size(), fpA, fpB);
188 }
189
191 << "Found connected faces " << mesh.faces()[f0I] << " and "
192 << mesh.faces()[f1I] << " sharing edge " << edgeI << endl
193 << "But none seems to be connected to face " << facei
194 << " vertices:" << f
195 << abort(FatalError);
196
197 return -1;
199 }
200}
201
202
203// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
204
205Foam::Ostream& Foam::operator<<
206(
207 Ostream& os,
208 const directionInfo& rhs
209)
210{
212 {
213 os << rhs.index_ << rhs.n_;
214 }
215 else
216 {
217 os.write
218 (
219 reinterpret_cast<const char*>(&rhs.index_),
220 sizeof(directionInfo)
221 );
223
224 os.check(FUNCTION_NAME);
225 return os;
226}
227
228
229Foam::Istream& Foam::operator>>
230(
231 Istream& is,
233)
234{
235 if (is.format() == IOstreamOption::ASCII)
236 {
237 is >> rhs.index_ >> rhs.n_;
238 }
239 else
240 {
241 // Packing of (label, vector) makes handling of non-native
242 // label/scalar non-trivial (#3412)
243 is.fatalCheckNativeSizes(FUNCTION_NAME);
244
245 is.read
246 (
247 reinterpret_cast<char*>(&rhs.index_),
248 sizeof(directionInfo)
249 );
250 }
251
252 is.check(FUNCTION_NAME);
253 return is;
254}
255
256// ************************************************************************* //
streamFormat format() const noexcept
Get the current stream format.
@ ASCII
"ascii" (normal default)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
virtual Ostream & write(const char c) override
Write character.
Definition OBJstream.C:69
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition UListI.H:454
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
Holds direction in which to split cell (in fact a local coordinate axes). Information is a label and ...
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.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
label & end() noexcept
The end (second/last) vertex label.
Definition edge.H:165
label start() const noexcept
The start (first) vertex label.
Definition edge.H:155
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
Cell-face mesh analysis engine.
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
label walkFace(const primitiveMesh &mesh, const label facei, const label startEdgeI, const label startVertI, const label nEdges)
Returns label of edge nEdges away from startEdge (in the direction.
Definition meshTools.C:596
void getEdgeFaces(const primitiveMesh &mesh, const label celli, const label edgeI, label &face0, label &face1)
Get faces on cell using edgeI. Throws error if no two found.
Definition meshTools.C:472
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
errorManip< error > abort(error &err)
Definition errorManip.H:139
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
labelList f(nPoints)
volScalarField & b
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299