boundaryPoints.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-2021 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 "boundaryPoints.H"
27 #include "polyMesh.H"
29 #include "treeBoundBox.H"
30 #include "treeDataFace.H"
31 #include "Time.H"
32 #include "meshTools.H"
33 #include "mappedPatchBase.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace sampledSets
40 {
41  defineTypeNameAndDebug(boundaryPoints, 0);
42  addToRunTimeSelectionTable(sampledSet, boundaryPoints, word);
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::sampledSets::boundaryPoints::calcSamples
50 (
51  DynamicList<point>& samplingPositions,
52  DynamicList<label>& samplingSegments,
53  DynamicList<label>& samplingCells,
54  DynamicList<label>& samplingFaces
55 ) const
56 {
57  // Get the patch IDs
58  const labelHashSet patchIDs(mesh().boundaryMesh().patchSet(patches_));
59 
60  // Construct a single list of all patch faces
61  label nPatchFaces = 0;
62  forAllConstIter(labelHashSet, patchIDs, iter)
63  {
64  const polyPatch& pp = mesh().boundaryMesh()[iter.key()];
65  nPatchFaces += pp.size();
66  }
67  labelList patchFaces(nPatchFaces);
68  nPatchFaces = 0;
69  forAllConstIter(labelHashSet, patchIDs, iter)
70  {
71  const polyPatch& pp = mesh().boundaryMesh()[iter.key()];
72  forAll(pp, i)
73  {
74  patchFaces[nPatchFaces++] = pp.start()+i;
75  }
76  }
77 
78  // Construct a processor-local bound box
79  treeBoundBox patchBB(point::max, point::min);
80  forAllConstIter(labelHashSet, patchIDs, iter)
81  {
82  const polyPatch& pp = mesh().boundaryMesh()[iter.key()];
83  const boundBox patchBb(pp.points(), pp.meshPoints(), false);
84  patchBB.min() = min(patchBB.min(), patchBb.min());
85  patchBB.max() = max(patchBB.max(), patchBb.max());
86  }
87  patchBB = patchBB.extend(1e-4);
88 
89  // Create a search tree
90  indexedOctree<treeDataFace> patchTree
91  (
92  treeDataFace // all information needed to search faces
93  (
94  false, // do not cache the bound box
95  mesh(),
96  patchFaces
97  ),
98  patchBB, // overall search box
99  8, // maximum number of levels
100  10, // how many elements per leaf
101  3.0 // in how many leaves is a shape on average
102  );
103 
104  // Force calculation of face-diagonal decomposition
105  (void)mesh().tetBasePtIs();
106 
107  // Generate the nearest patch information for each sampling point
108  List<mappedPatchBase::nearInfo> nearest(points_.size());
109  forAll(points_, sampleI)
110  {
111  const point& sample = points_[sampleI];
112 
113  pointIndexHit& nearHit = nearest[sampleI].first();
114  scalar& nearDist = nearest[sampleI].second().first();
115  label& nearProc = nearest[sampleI].second().second();
116 
117  // Find the nearest
118  if (patchFaces.size())
119  {
120  nearHit = patchTree.findNearest(sample, sqr(maxDistance_));
121  }
122  else
123  {
124  nearHit.setMiss();
125  }
126 
127  // Fill in the information
128  if (nearHit.hit())
129  {
130  nearHit.setIndex(patchFaces[nearHit.index()]);
131  nearDist = magSqr(nearHit.hitPoint() - sample);
132  nearProc = Pstream::myProcNo();
133  }
134  else
135  {
136  nearHit.setIndex(-1);
137  nearDist = Foam::sqr(great);
138  nearProc = Pstream::myProcNo();
139  }
140  }
141 
142  // Reduce to get the nearest patch locations globally
143  Pstream::listCombineGather(nearest, mappedPatchBase::nearestEqOp());
145 
146  // Dump connecting lines from the sampling points to the hit locations
147  if (debug && Pstream::master())
148  {
149  OFstream str(mesh().time().path() / name() + "_nearest.obj");
150 
151  label verti = 0;
152 
153  forAll(nearest, i)
154  {
155  if (nearest[i].first().hit())
156  {
157  meshTools::writeOBJ(str, points_[i]);
158  verti++;
159  meshTools::writeOBJ(str, nearest[i].first().hitPoint());
160  verti++;
161  str << "l " << verti - 1 << ' ' << verti << nl;
162  }
163  }
164  }
165 
166  // Store the sampling locations on the nearest processor
167  forAll(nearest, sampleI)
168  {
169  const pointIndexHit& nearHit = nearest[sampleI].first();
170  const label& nearProc = nearest[sampleI].second().second();
171 
172  if (nearHit.hit())
173  {
174  if (nearProc == Pstream::myProcNo())
175  {
176  label facei = nearHit.index();
177 
178  samplingPositions.append(nearHit.hitPoint());
179  samplingSegments.append(sampleI);
180  samplingCells.append(mesh().faceOwner()[facei]);
181  samplingFaces.append(facei);
182  }
183  }
184  else
185  {
187  << "Unable to find location on patches " << patches_
188  << " for the point " << points_[sampleI]
189  << " within a distance of " << maxDistance_ << endl;
190  }
191  }
192 }
193 
194 
195 void Foam::sampledSets::boundaryPoints::genSamples()
196 {
197  DynamicList<point> samplingPositions;
198  DynamicList<label> samplingSegments;
199  DynamicList<label> samplingCells;
200  DynamicList<label> samplingFaces;
201 
202  calcSamples
203  (
204  samplingPositions,
205  samplingSegments,
206  samplingCells,
207  samplingFaces
208  );
209 
210  samplingPositions.shrink();
211  samplingSegments.shrink();
212  samplingCells.shrink();
213  samplingFaces.shrink();
214 
215  setSamples
216  (
217  samplingPositions,
218  samplingSegments,
219  samplingCells,
220  samplingFaces
221  );
222 }
223 
224 
225 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
226 
228 (
229  const word& name,
230  const polyMesh& mesh,
231  const meshSearch& searchEngine,
232  const dictionary& dict
233 )
234 :
235  sampledSet(name, mesh, searchEngine, dict),
236  points_(dict.lookup("points")),
237  patches_(dict.lookup("patches")),
238  maxDistance_(dict.lookup<scalar>("maxDistance"))
239 {
240  genSamples();
241 }
242 
243 
244 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
245 
247 {}
248 
249 
250 // ************************************************************************* //
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:57
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:38
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
boundaryPoints(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
fvMesh & mesh
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Macros for easy insertion into run-time selection tables.
static const layerAndWeight max
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:61
A class for handling words, derived from string.
Definition: word.H:59
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
List< label > labelList
A List of labels.
Definition: labelList.H:56
addToRunTimeSelectionTable(sampledSet, arcUniform, word)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual ~boundaryPoints()
Destructor.
static const char nl
Definition: Ostream.H:260
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
vector point
Point is a vector.
Definition: point.H:41
#define WarningInFunction
Report a warning using Foam::Warning.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
defineTypeNameAndDebug(arcUniform, 0)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864