PatchToolsCheck.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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "PatchTools.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 template<class FaceList, class PointField>
34 bool
36 (
38  const bool report,
39  labelHashSet* setPtr
40 )
41 {
42  typedef typename PrimitivePatch<FaceList, PointField>::face_type FaceType;
43 
44  bool foundError = false;
45 
46  // Check edge normals, face normals, point normals.
47  forAll(p.faceEdges(), facei)
48  {
49  const labelList& edgeLabels = p.faceEdges()[facei];
50  bool valid = true;
51 
52  if (edgeLabels.size() < 3)
53  {
54  if (report)
55  {
56  Info<< "Face[" << facei << "] " << p[facei]
57  << " has fewer than 3 edges. Edges: " << edgeLabels
58  << endl;
59  }
60  valid = false;
61  }
62  else
63  {
64  forAll(edgeLabels, i)
65  {
66  if (edgeLabels[i] < 0 || edgeLabels[i] >= p.nEdges())
67  {
68  if (report)
69  {
70  Info<< "edge number " << edgeLabels[i]
71  << " on face " << facei
72  << " out-of-range\n"
73  << "This usually means the input surface has "
74  << "edges with more than 2 faces connected."
75  << endl;
76  }
77  valid = false;
78  }
79  }
80  }
81 
82  if (!valid)
83  {
84  foundError = true;
85  continue;
86  }
87 
88 
89  //
90  //- Compute normal from 3 points, use the first as the origin
91  // minor warpage should not be a problem
92  const FaceType& f = p[facei];
93  const point& p0 = p.points()[f[0]];
94  const point& p1 = p.points()[f[1]];
95  const point& p2 = p.points()[f.last()];
96 
97  const vector pointNormal((p1 - p0) ^ (p2 - p0));
98  if ((pointNormal & p.faceNormals()[facei]) < 0)
99  {
100  foundError = false;
101 
102  if (report)
103  {
104  Info
105  << "Normal calculated from points inconsistent"
106  << " with faceNormal" << nl
107  << "face: " << f << nl
108  << "points: " << p0 << ' ' << p1 << ' ' << p2 << nl
109  << "pointNormal:" << pointNormal << nl
110  << "faceNormal:" << p.faceNormals()[facei] << endl;
111  }
112  }
113  }
114 
115 
116  forAll(p.edges(), edgeI)
117  {
118  const edge& e = p.edges()[edgeI];
119  const labelList& neighbouringFaces = p.edgeFaces()[edgeI];
120 
121  if (neighbouringFaces.size() == 2)
122  {
123  // we use localFaces() since edges() are LOCAL
124  // these are both already available
125  const FaceType& faceA = p.localFaces()[neighbouringFaces[0]];
126  const FaceType& faceB = p.localFaces()[neighbouringFaces[1]];
127 
128  // If the faces are correctly oriented, the edges must go in
129  // different directions on connected faces.
130  if (faceA.edgeDirection(e) == faceB.edgeDirection(e))
131  {
132  if (report)
133  {
134  Info<< "face orientation incorrect." << nl
135  << "localEdge[" << edgeI << "] " << e
136  << " between faces:" << nl
137  << " face[" << neighbouringFaces[0] << "] "
138  << p[neighbouringFaces[0]]
139  << " localFace: " << faceA
140  << nl
141  << " face[" << neighbouringFaces[1] << "] "
142  << p[neighbouringFaces[1]]
143  << " localFace: " << faceB
144  << endl;
145  }
146  if (setPtr)
147  {
148  setPtr->insert(edgeI);
149  }
150 
151  foundError = true;
152  }
153  }
154  else if (neighbouringFaces.size() != 1)
155  {
156  if (report)
157  {
158  Info
159  << "Wrong number of edge neighbours." << nl
160  << "edge[" << edgeI << "] " << e
161  << " with points:" << p.localPoints()[e.start()]
162  << ' ' << p.localPoints()[e.end()]
163  << " has neighbouringFaces:" << neighbouringFaces << endl;
164  }
165 
166  if (setPtr)
167  {
168  setPtr->insert(edgeI);
169  }
170 
171  foundError = true;
172  }
173  }
174 
175  return foundError;
176 }
177 
178 
179 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
PatchTools.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::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::HashSet< label, Hash< label > >
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::nl
constexpr char nl
Definition: Ostream.H:404
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
Foam::PrimitivePatch::face_type
std::remove_reference< FaceList >::type::value_type face_type
The face type.
Definition: PrimitivePatch.H:90
Foam::PatchTools::checkOrientation
static bool checkOrientation(const PrimitivePatch< FaceList, PointField > &, const bool report=false, labelHashSet *marked=0)
Check for orientation issues.
Definition: PatchToolsCheck.C:36
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79