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-2022 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 "RemoteData.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace sampledSets
40 {
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<RemoteData<Tuple2<scalar, point>>> nearest(points_.size());
109  forAll(points_, sampleI)
110  {
111  const point& sample = points_[sampleI];
112 
113  const pointIndexHit pih =
114  patchFaces.size()
115  ? patchTree.findNearest(sample, sqr(maxDistance_))
116  : pointIndexHit();
117 
118  if (pih.hit())
119  {
120  nearest[sampleI].proci = Pstream::myProcNo();
121  nearest[sampleI].elementi = patchFaces[pih.index()];
122  nearest[sampleI].data.first() = magSqr(pih.hitPoint() - sample);
123  nearest[sampleI].data.second() = pih.hitPoint();
124  }
125  }
126 
127  // Reduce to get the nearest patch locations globally
129  (
130  nearest,
131  RemoteData<Tuple2<scalar, point>>::smallestFirstEqOp()
132  );
134 
135  // Dump connecting lines from the sampling points to the hit locations
136  if (debug && Pstream::master())
137  {
138  OFstream str(mesh().time().path() / name() + "_nearest.obj");
139 
140  label verti = 0;
141 
142  forAll(nearest, sampleI)
143  {
144  if (nearest[sampleI].proci != -1)
145  {
146  meshTools::writeOBJ(str, points_[sampleI]);
147  verti++;
148  meshTools::writeOBJ(str, nearest[sampleI].data.second());
149  verti++;
150  str << "l " << verti - 1 << ' ' << verti << nl;
151  }
152  }
153  }
154 
155  // Store the sampling locations on the nearest processor
156  forAll(nearest, sampleI)
157  {
158  if (nearest[sampleI].proci != -1)
159  {
160  if (nearest[sampleI].proci == Pstream::myProcNo())
161  {
162  const label facei = nearest[sampleI].elementi;
163 
164  samplingPositions.append(nearest[sampleI].data.second());
165  samplingSegments.append(sampleI);
166  samplingCells.append(mesh().faceOwner()[facei]);
167  samplingFaces.append(facei);
168  }
169  }
170  else
171  {
173  << "Unable to find location on patches " << patches_
174  << " for the point " << points_[sampleI]
175  << " within a distance of " << maxDistance_ << endl;
176  }
177  }
178 }
179 
180 
181 void Foam::sampledSets::boundaryPoints::genSamples()
182 {
183  DynamicList<point> samplingPositions;
184  DynamicList<label> samplingSegments;
185  DynamicList<label> samplingCells;
186  DynamicList<label> samplingFaces;
187 
188  calcSamples
189  (
190  samplingPositions,
191  samplingSegments,
192  samplingCells,
193  samplingFaces
194  );
195 
196  samplingPositions.shrink();
197  samplingSegments.shrink();
198  samplingCells.shrink();
199  samplingFaces.shrink();
200 
201  setSamples
202  (
203  samplingPositions,
204  samplingSegments,
205  samplingCells,
206  samplingFaces
207  );
208 }
209 
210 
211 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
212 
214 (
215  const word& name,
216  const polyMesh& mesh,
217  const meshSearch& searchEngine,
218  const dictionary& dict
219 )
220 :
221  sampledSet(name, mesh, searchEngine, dict),
222  points_(dict.lookup("points")),
223  patches_(dict.lookup("patches")),
224  maxDistance_(dict.lookup<scalar>("maxDistance"))
225 {
226  genSamples();
227 }
228 
229 
230 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
231 
233 {}
234 
235 
236 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static const Form max
Definition: VectorSpace.H:115
static const Form min
Definition: VectorSpace.H:116
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:58
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:1078
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:64
const word & name() const
Access the name.
Definition: sampledSet.H:204
const polyMesh & mesh() const
Access the mesh.
Definition: sampledSet.H:210
Specified point samples within patches.
boundaryPoints(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
virtual ~boundaryPoints()
Destructor.
Set of sets to sample. Call sampledSets.write() to sample&write files.
A class for handling words, derived from string.
Definition: word.H:62
#define WarningInFunction
Report a warning using Foam::Warning.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
defineTypeNameAndDebug(arcUniform, 0)
addToRunTimeSelectionTable(sampledSet, arcUniform, word)
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:105
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
vector point
Point is a vector.
Definition: point.H:41
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict