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-2020 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>& samplingPts,
52  DynamicList<label>& samplingCells,
53  DynamicList<label>& samplingFaces,
54  DynamicList<label>& samplingSegments,
55  DynamicList<scalar>& samplingCurveDist
56 ) const
57 {
58  // Get the patch IDs
59  const labelHashSet patchIDs(mesh().boundaryMesh().patchSet(patches_));
60 
61  // Construct a single list of all patch faces
62  label nPatchFaces = 0;
63  forAllConstIter(labelHashSet, patchIDs, iter)
64  {
65  const polyPatch& pp = mesh().boundaryMesh()[iter.key()];
66  nPatchFaces += pp.size();
67  }
68  labelList patchFaces(nPatchFaces);
69  nPatchFaces = 0;
70  forAllConstIter(labelHashSet, patchIDs, iter)
71  {
72  const polyPatch& pp = mesh().boundaryMesh()[iter.key()];
73  forAll(pp, i)
74  {
75  patchFaces[nPatchFaces++] = pp.start()+i;
76  }
77  }
78 
79  // Construct a processor-local bound box
80  treeBoundBox patchBB(point::max, point::min);
81  forAllConstIter(labelHashSet, patchIDs, iter)
82  {
83  const polyPatch& pp = mesh().boundaryMesh()[iter.key()];
84  const boundBox patchBb(pp.points(), pp.meshPoints(), false);
85  patchBB.min() = min(patchBB.min(), patchBb.min());
86  patchBB.max() = max(patchBB.max(), patchBb.max());
87  }
88  patchBB = patchBB.extend(1e-4);
89 
90  // Create a search tree
91  indexedOctree<treeDataFace> patchTree
92  (
93  treeDataFace // all information needed to search faces
94  (
95  false, // do not cache the bound box
96  mesh(),
97  patchFaces
98  ),
99  patchBB, // overall search box
100  8, // maximum number of levels
101  10, // how many elements per leaf
102  3.0 // in how many leaves is a shape on average
103  );
104 
105  // Force calculation of face-diagonal decomposition
106  (void)mesh().tetBasePtIs();
107 
108  // Generate the nearest patch information for each sampling point
109  List<mappedPatchBase::nearInfo> nearest(points_.size());
110  forAll(points_, sampleI)
111  {
112  const point& sample = points_[sampleI];
113 
114  pointIndexHit& nearHit = nearest[sampleI].first();
115  scalar& nearDist = nearest[sampleI].second().first();
116  label& nearProc = nearest[sampleI].second().second();
117 
118  // Find the nearest
119  if (patchFaces.size())
120  {
121  nearHit = patchTree.findNearest(sample, sqr(maxDistance_));
122  }
123  else
124  {
125  nearHit.setMiss();
126  }
127 
128  // Fill in the information
129  if (nearHit.hit())
130  {
131  nearHit.setIndex(patchFaces[nearHit.index()]);
132  nearDist = magSqr(nearHit.hitPoint() - sample);
133  nearProc = Pstream::myProcNo();
134  }
135  else
136  {
137  nearHit.setIndex(-1);
138  nearDist = Foam::sqr(great);
139  nearProc = Pstream::myProcNo();
140  }
141  }
142 
143  // Reduce to get the nearest patch locations globally
144  Pstream::listCombineGather(nearest, mappedPatchBase::nearestEqOp());
146 
147  // Dump connecting lines from the sampling points to the hit locations
148  if (debug && Pstream::master())
149  {
150  OFstream str(mesh().time().path() / name() + "_nearest.obj");
151 
152  label verti = 0;
153 
154  forAll(nearest, i)
155  {
156  if (nearest[i].first().hit())
157  {
158  meshTools::writeOBJ(str, points_[i]);
159  verti++;
160  meshTools::writeOBJ(str, nearest[i].first().hitPoint());
161  verti++;
162  str << "l " << verti - 1 << ' ' << verti << nl;
163  }
164  }
165  }
166 
167  // Store the sampling locations on the nearest processor
168  forAll(nearest, sampleI)
169  {
170  const pointIndexHit& nearHit = nearest[sampleI].first();
171  const label& nearProc = nearest[sampleI].second().second();
172 
173  if (nearHit.hit())
174  {
175  if (nearProc == Pstream::myProcNo())
176  {
177  label facei = nearHit.index();
178 
179  samplingPts.append(nearHit.hitPoint());
180  samplingCells.append(mesh().faceOwner()[facei]);
181  samplingFaces.append(facei);
182  samplingSegments.append(0);
183  samplingCurveDist.append(sampleI);
184  }
185  }
186  else
187  {
189  << "Unable to find location on patches " << patches_
190  << " for the point " << points_[sampleI]
191  << " within a distance of " << maxDistance_ << endl;
192  }
193  }
194 }
195 
196 
197 void Foam::sampledSets::boundaryPoints::genSamples()
198 {
199  DynamicList<point> samplingPts;
200  DynamicList<label> samplingCells;
201  DynamicList<label> samplingFaces;
202  DynamicList<label> samplingSegments;
203  DynamicList<scalar> samplingCurveDist;
204 
205  calcSamples
206  (
207  samplingPts,
208  samplingCells,
209  samplingFaces,
210  samplingSegments,
211  samplingCurveDist
212  );
213 
214  samplingPts.shrink();
215  samplingCells.shrink();
216  samplingFaces.shrink();
217  samplingSegments.shrink();
218  samplingCurveDist.shrink();
219 
220  setSamples
221  (
222  samplingPts,
223  samplingCells,
224  samplingFaces,
225  samplingSegments,
226  samplingCurveDist
227  );
228 }
229 
230 
231 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
232 
234 (
235  const word& name,
236  const polyMesh& mesh,
237  const meshSearch& searchEngine,
238  const dictionary& dict
239 )
240 :
241  sampledSet(name, mesh, searchEngine, dict),
242  points_(dict.lookup("points")),
243  patches_(dict.lookup("patches")),
244  maxDistance_(dict.lookup<scalar>("maxDistance"))
245 {
246  genSamples();
247 
248  if (debug)
249  {
250  write(Info);
251  }
252 }
253 
254 
255 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
256 
258 {}
259 
260 
261 // ************************************************************************* //
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:57
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
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
boundaryPoints(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
const word & name() const
Return the name of this functionObject.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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
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.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
dynamicFvMesh & mesh
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
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
addToRunTimeSelectionTable(sampledSet, arcUniform, word)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual ~boundaryPoints()
Destructor.
static const char nl
Definition: Ostream.H:260
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
vector point
Point is a vector.
Definition: point.H:41
#define WarningInFunction
Report a warning using Foam::Warning.
virtual bool write()
Sample and write.
Definition: sampledSets.C:228
messageStream Info
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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:844
fileName path(UMean.rootPath()/UMean.caseName()/functionObjects::writeFile::outputPrefix/"graphs"/UMean.instance())
T & first()
Return reference to the first element of the list.
Definition: UPtrListI.H:43