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