Loading...
Searching...
No Matches
cellDistFuncsTemplates.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,2024 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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 "cellDistFuncs.H"
29#include "polyMesh.H"
30#include "polyBoundaryMesh.H"
31
32
33// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34
35template<class Type>
37{
39
41
42 forAll(bMesh, patchi)
43 {
44 if (isA<Type>(bMesh[patchi]))
45 {
46 patchIDs.insert(patchi);
47 }
48 }
49 return patchIDs;
50}
51
52
53template<class PatchType>
55(
56 const point& p,
57 const PatchType& patch,
58 const labelUList& wallFaces,
59 label& minFacei
60) const
61{
62 // Return smallest true distance from p to any of wallFaces.
63 // Note that even if normal hits face we still check other faces.
64
65 const pointField& points = patch.points();
66
67 scalar minDist = GREAT;
68 minFacei = -1;
69
70 for (const label patchFacei : wallFaces)
71 {
72 const pointHit curHit = patch[patchFacei].nearestPoint(p, points);
73
74 if (curHit.distance() < minDist)
75 {
76 minDist = curHit.distance();
77 minFacei = patchFacei;
78 }
79 }
80
81 return minDist;
82}
83
84
85template<class PatchType>
87(
88 const PatchType& patch,
89 const label patchFacei,
90 DynamicList<label>& neighbours
91) const
92{
93 // Get point neighbours of facei (including facei). Returns number of faces.
94 // Note: does not allocate storage but does use linear search to determine
95 // uniqueness. For polygonal faces this might be quite inefficient.
96
97 neighbours.clear();
98
99 // Add myself
100 neighbours.append(patchFacei);
101
102 // Add all face neighbours
103 const labelList& faceNeighbours = patch.faceFaces()[patchFacei];
104
105 for (const label nbr : faceNeighbours)
106 {
107 neighbours.push_uniq(nbr);
108 }
109
110 // Add all point-only neighbours by linear searching in edge neighbours.
111 // Assumes that point-only neighbours are not using multiple points on
112 // face.
113
114 const face& f = patch.localFaces()[patchFacei];
115
116 forAll(f, fp)
117 {
118 label pointi = f[fp];
119
120 const labelList& pointNbs = patch.pointFaces()[pointi];
121
122 for (const label facei : pointNbs)
123 {
124 // Check for facei in edge-neighbours part of neighbours
125 neighbours.push_uniq(facei);
126 }
127 }
128
129
130 if (debug)
131 {
132 // Check for duplicates
133
134 // Use hashSet to determine nbs.
135 labelHashSet nbs(4*f.size());
136
137 forAll(f, fp)
138 {
139 const labelList& pointNbs = patch.pointFaces()[f[fp]];
140 nbs.insert(pointNbs);
141 }
142
143 // Subtract ours.
144 for (const label nb : neighbours)
145 {
146 if (!nbs.found(nb))
147 {
149 << "getPointNeighbours : patchFacei:" << patchFacei
150 << " verts:" << f << endl;
151
152 forAll(f, fp)
153 {
155 << "point:" << f[fp] << " pointFaces:"
156 << patch.pointFaces()[f[fp]] << endl;
157 }
158
159 for (const label facei : neighbours)
160 {
162 << "fast nbr:" << facei
163 << endl;
164 }
165
167 << "Problem: fast pointNeighbours routine included " << nb
168 << " which is not in proper neighbour list " << nbs.toc()
169 << abort(FatalError);
170 }
171 nbs.erase(nb);
172 }
173
174 if (nbs.size())
175 {
177 << "Problem: fast pointNeighbours routine did not find "
178 << nbs.toc() << abort(FatalError);
179 }
180 }
181}
182
183
184// ************************************************************************* //
labelList patchIDs
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void append(const T &val)
Copy append an element to the end of this list.
label push_uniq(const T &val)
Append an element if not already in the list.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition HashTable.C:141
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition HashTable.C:489
scalar distance() const noexcept
Return distance to hit.
Definition pointHit.H:169
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
labelHashSet getPatchIDs() const
Get patchIDs of/derived off certain type (e.g. 'processorPolyPatch').
const polyMesh & mesh() const
Access mesh.
scalar smallestDist(const point &p, const PatchType &patch, const labelUList &wallFaces, label &patchFacei) const
Calculate smallest true distance (and patch face index).
void getPointNeighbours(const PatchType &, const label patchFacei, DynamicList< label > &) const
Get faces sharing point with face on patch.
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
errorManip< error > abort(error &err)
Definition errorManip.H:139
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points).
Definition bMesh.H:39
vectorField pointField
pointField is a vectorField.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299