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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "directionInfo.H"
29 #include "polyMesh.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 // Find edge among edgeLabels that uses v0 and v1
34 Foam::label Foam::directionInfo::findEdge
35 (
36  const primitiveMesh& mesh,
37  const labelList& edgeLabels,
38  const label v1,
39  const label v0
40 )
41 {
42  forAll(edgeLabels, edgeLabelI)
43  {
44  label edgeI = edgeLabels[edgeLabelI];
45 
46  if (mesh.edges()[edgeI] == edge(v0, v1))
47  {
48  return edgeI;
49  }
50  }
51 
53  << "Cannot find an edge among " << edgeLabels << endl
54  << "that uses vertices " << v0
55  << " and " << v1
56  << abort(FatalError);
57 
58  return -1;
59 }
60 
61 
62 Foam::label Foam::directionInfo::lowest
63 (
64  const label size,
65  const label a,
66  const label b
67 )
68 {
69  // Get next point
70  label a1 = (a + 1) % size;
71 
72  if (a1 == b)
73  {
74  return a;
75  }
76  else
77  {
78  label b1 = (b + 1) % size;
79 
80  if (b1 != a)
81  {
83  << "Problem : a:" << a << " b:" << b << " size:" << size
84  << abort(FatalError);
85  }
86 
87  return b;
88  }
89 }
90 
91 
92 // Have edge on hex cell. Find corresponding edge on face. Return -1 if none
93 // found.
95 (
96  const primitiveMesh& mesh,
97  const label celli,
98  const label facei,
99  const label edgeI
100 )
101 {
102  if ((edgeI < 0) || (edgeI >= mesh.nEdges()))
103  {
105  << "Illegal edge label:" << edgeI
106  << " when projecting cut edge from cell " << celli
107  << " to face " << facei
108  << abort(FatalError);
109  }
110 
111  const edge& e = mesh.edges()[edgeI];
112 
113  const face& f = mesh.faces()[facei];
114 
115  // edgeI is either
116  // - in facei. Convert into index in face.
117  // - connected (but not in) to face. Return -1.
118  // - in face opposite facei. Convert into index in face.
119 
120  label fpA = f.find(e.start());
121  label fpB = f.find(e.end());
122 
123  if (fpA != -1)
124  {
125  if (fpB != -1)
126  {
127  return lowest(f.size(), fpA, fpB);
128  }
129  else
130  {
131  // e.start() in face, e.end() not
132  return -1;
133  }
134  }
135  else
136  {
137  if (fpB != -1)
138  {
139  // e.end() in face, e.start() not
140  return -1;
141  }
142  else
143  {
144  // Both not in face.
145  // e is on opposite face. Determine corresponding edge on this face:
146  // - determine two faces using edge (one is the opposite face,
147  // one is 'side' face
148  // - walk on both these faces to opposite edge
149  // - check if this opposite edge is on facei
150 
151  label f0I, f1I;
152 
153  meshTools::getEdgeFaces(mesh, celli, edgeI, f0I, f1I);
154 
155  // Walk to opposite edge on face f0
156  label edge0I =
157  meshTools::walkFace(mesh, f0I, edgeI, e.start(), 2);
158 
159  // Check if edge on facei.
160 
161  const edge& e0 = mesh.edges()[edge0I];
162 
163  fpA = f.find(e0.start());
164  fpB = f.find(e0.end());
165 
166  if ((fpA != -1) && (fpB != -1))
167  {
168  return lowest(f.size(), fpA, fpB);
169  }
170 
171  // Face0 is doesn't have an edge on facei (so must be the opposite
172  // face) so try face1.
173 
174  // Walk to opposite edge on face f1
175  label edge1I =
176  meshTools::walkFace(mesh, f1I, edgeI, e.start(), 2);
177 
178  // Check if edge on facei.
179  const edge& e1 = mesh.edges()[edge1I];
180 
181  fpA = f.find(e1.start());
182  fpB = f.find(e1.end());
183 
184  if ((fpA != -1) && (fpB != -1))
185  {
186  return lowest(f.size(), fpA, fpB);
187  }
188 
190  << "Found connected faces " << mesh.faces()[f0I] << " and "
191  << mesh.faces()[f1I] << " sharing edge " << edgeI << endl
192  << "But none seems to be connected to face " << facei
193  << " vertices:" << f
194  << abort(FatalError);
195 
196  return -1;
197  }
198  }
199 }
200 
201 
202 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
203 
204 Foam::Ostream& Foam::operator<<
205 (
206  Foam::Ostream& os,
207  const Foam::directionInfo& wDist
208 )
209 {
210  if (os.format() == IOstream::ASCII)
211  {
212  os << wDist.index_ << wDist.n_;
213  }
214  else
215  {
216  os.write
217  (
218  reinterpret_cast<const char*>(&wDist.index_),
219  sizeof(directionInfo)
220  );
221  }
222 
223  os.check(FUNCTION_NAME);
224  return os;
225 }
226 
227 
229 {
230  if (is.format() == IOstream::ASCII)
231  {
232  is >> wDist.index_ >> wDist.n_;
233  }
234  else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
235  {
236  // Non-native label or scalar size
237  is.beginRawRead();
238 
239  readRawLabel(is, &wDist.index_);
240  readRawScalar(is, wDist.n_.data(), vector::nComponents);
241 
242  is.endRawRead();
243  }
244  else
245  {
246  is.read
247  (
248  reinterpret_cast<char*>(&wDist.index_),
249  sizeof(directionInfo)
250  );
251  }
252 
253  is.check(FUNCTION_NAME);
254  return is;
255 }
256 
257 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::directionInfo::edgeToFaceIndex
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:95
directionInfo.H
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::IOstreamOption::format
streamFormat format() const noexcept
Get the current stream format.
Definition: IOstreamOption.H:273
Foam::primitiveMesh::nEdges
label nEdges() const
Number of mesh edges.
Definition: primitiveMeshI.H:67
Foam::VectorSpace::data
Cmpt * data() noexcept
Return pointer to the first data element.
Definition: VectorSpaceI.H:198
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:228
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
polyMesh.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::meshTools::walkFace
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:603
Foam::label
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:62
Foam::Istream::endRawRead
virtual bool endRawRead()=0
End of low-level raw binary read.
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::IOstream::checkLabelSize
std::enable_if< std::is_integral< T >::value, bool >::type checkLabelSize() const
Definition: IOstream.H:283
Foam::meshTools::getEdgeFaces
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:479
Foam::directionInfo
Holds direction in which to split cell (in fact a local coordinate axes). Information is a label and ...
Definition: directionInfo.H:85
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:51
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::IOstreamOption::ASCII
"ascii"
Definition: IOstreamOption.H:66
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::Istream::beginRawRead
virtual bool beginRawRead()=0
Start of low-level raw binary read.
Foam::IOstream::checkScalarSize
std::enable_if< std::is_floating_point< T >::value, bool >::type checkScalarSize() const
Definition: IOstream.H:292
f
labelList f(nPoints)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:261
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
Foam::readRawLabel
label readRawLabel(Istream &is)
Read raw label from binary stream.
Definition: label.C:46
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::nComponents
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:101
Foam::Istream::read
virtual Istream & read(token &)=0
Return next token from stream.
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:78