directionInfo.C
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) 2011-2018 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 "directionInfo.H"
27 #include "polyMesh.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 // Find edge among edgeLabels that uses v0 and v1
32 Foam::label Foam::directionInfo::findEdge
33 (
34  const primitiveMesh& mesh,
35  const labelList& edgeLabels,
36  const label v1,
37  const label v0
38 )
39 {
40  forAll(edgeLabels, edgeLabelI)
41  {
42  label edgeI = edgeLabels[edgeLabelI];
43 
44  if (mesh.edges()[edgeI] == edge(v0, v1))
45  {
46  return edgeI;
47  }
48  }
49 
51  << "Cannot find an edge among " << edgeLabels << endl
52  << "that uses vertices " << v0
53  << " and " << v1
54  << abort(FatalError);
55 
56  return -1;
57 }
58 
59 
60 Foam::label Foam::directionInfo::lowest
61 (
62  const label size,
63  const label a,
64  const label b
65 )
66 {
67  // Get next point
68  label a1 = (a + 1) % size;
69 
70  if (a1 == b)
71  {
72  return a;
73  }
74  else
75  {
76  label b1 = (b + 1) % size;
77 
78  if (b1 != a)
79  {
81  << "Problem : a:" << a << " b:" << b << " size:" << size
82  << abort(FatalError);
83  }
84 
85  return b;
86  }
87 }
88 
89 
90 // Have edge on hex cell. Find corresponding edge on face. Return -1 if none
91 // found.
93 (
94  const primitiveMesh& mesh,
95  const label celli,
96  const label facei,
97  const label edgeI
98 )
99 {
100  if ((edgeI < 0) || (edgeI >= mesh.nEdges()))
101  {
103  << "Illegal edge label:" << edgeI
104  << " when projecting cut edge from cell " << celli
105  << " to face " << facei
106  << abort(FatalError);
107  }
108 
109  const edge& e = mesh.edges()[edgeI];
110 
111  const face& f = mesh.faces()[facei];
112 
113  // edgeI is either
114  // - in facei. Convert into index in face.
115  // - connected (but not in) to face. Return -1.
116  // - in face opposite facei. Convert into index in face.
117 
118  label fpA = findIndex(f, e.start());
119  label fpB = findIndex(f, e.end());
120 
121  if (fpA != -1)
122  {
123  if (fpB != -1)
124  {
125  return lowest(f.size(), fpA, fpB);
126  }
127  else
128  {
129  // e.start() in face, e.end() not
130  return -1;
131  }
132  }
133  else
134  {
135  if (fpB != -1)
136  {
137  // e.end() in face, e.start() not
138  return -1;
139  }
140  else
141  {
142  // Both not in face.
143  // e is on opposite face. Determine corresponding edge on this face:
144  // - determine two faces using edge (one is the opposite face,
145  // one is 'side' face
146  // - walk on both these faces to opposite edge
147  // - check if this opposite edge is on facei
148 
149  label f0I, f1I;
150 
151  meshTools::getEdgeFaces(mesh, celli, edgeI, f0I, f1I);
152 
153  // Walk to opposite edge on face f0
154  label edge0I =
155  meshTools::walkFace(mesh, f0I, edgeI, e.start(), 2);
156 
157  // Check if edge on facei.
158 
159  const edge& e0 = mesh.edges()[edge0I];
160 
161  fpA = findIndex(f, e0.start());
162  fpB = findIndex(f, e0.end());
163 
164  if ((fpA != -1) && (fpB != -1))
165  {
166  return lowest(f.size(), fpA, fpB);
167  }
168 
169  // Face0 is doesn't have an edge on facei (so must be the opposite
170  // face) so try face1.
171 
172  // Walk to opposite edge on face f1
173  label edge1I =
174  meshTools::walkFace(mesh, f1I, edgeI, e.start(), 2);
175 
176  // Check if edge on facei.
177  const edge& e1 = mesh.edges()[edge1I];
178 
179  fpA = findIndex(f, e1.start());
180  fpB = findIndex(f, e1.end());
181 
182  if ((fpA != -1) && (fpB != -1))
183  {
184  return lowest(f.size(), fpA, fpB);
185  }
186 
188  << "Found connected faces " << mesh.faces()[f0I] << " and "
189  << mesh.faces()[f1I] << " sharing edge " << edgeI << endl
190  << "But none seems to be connected to face " << facei
191  << " vertices:" << f
192  << abort(FatalError);
193 
194  return -1;
195  }
196  }
197 }
198 
199 
200 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
201 
202 Foam::Ostream& Foam::operator<<
203 (
204  Foam::Ostream& os,
205  const Foam::directionInfo& wDist
206 )
207 {
208  if (os.format() == IOstream::ASCII)
209  {
210  os << wDist.index_ << wDist.n_;
211  }
212  else
213  {
214  os.write
215  (
216  reinterpret_cast<const char*>(&wDist.index_),
217  sizeof(directionInfo)
218  );
219  }
220 
221  // Check state of Ostream
222  os.check("Ostream& operator<<(Ostream&, const directionInfo&)");
223  return os;
224 
225 }
226 
227 
229 {
230  if (is.format() == IOstream::ASCII)
231  {
232  is >> wDist.index_ >> wDist.n_;
233  }
234  else
235  {
236  is.read
237  (
238  reinterpret_cast<char*>(&wDist.index_),
239  sizeof(directionInfo)
240  );
241  }
242 
243  // Check state of Istream
244  is.check("Istream& operator>>(Istream&, directionInfo&)");
245  return is;
246 }
247 
248 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Holds information (coordinate and normal) regarding nearest wall point.
iterator end()
Return an iterator to end traversing the UList.
Definition: UListI.H:224
virtual Istream & read(token &)=0
Return next token from stream.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Istream & operator>>(Istream &, directionInfo &)
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:93
List< label > labelList
A List of labels.
Definition: labelList.H:56
streamFormat format() const
Return current stream format.
Definition: IOstream.H:374
errorManip< error > abort(error &err)
Definition: errorManip.H:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
labelList f(nPoints)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
Holds direction in which to split cell (in fact a local coordinate axes). Information is a label and ...
Definition: directionInfo.H:85
label nEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
label end() const
Return end vertex label.
Definition: edgeI.H:92
virtual const faceList & faces() const =0
Return faces.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
label walkFace(const primitiveMesh &, 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:587
void getEdgeFaces(const primitiveMesh &, 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:458
label start() const
Return start vertex label.
Definition: edgeI.H:81