patchCloudSet.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  Copyright (C) 2017-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "patchCloudSet.H"
30 #include "polyMesh.H"
32 #include "treeBoundBox.H"
33 #include "treeDataFace.H"
34 #include "Time.H"
35 #include "meshTools.H"
36 // For 'nearInfo' helper class only
37 #include "mappedPatchBase.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(patchCloudSet, 0);
44  addToRunTimeSelectionTable(sampledSet, patchCloudSet, word);
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::patchCloudSet::calcSamples
51 (
52  DynamicList<point>& samplingPts,
53  DynamicList<label>& samplingCells,
54  DynamicList<label>& samplingFaces,
55  DynamicList<label>& samplingSegments,
56  DynamicList<scalar>& samplingCurveDist
57 ) const
58 {
59  DebugInfo << "patchCloudSet : sampling on patches :" << endl;
60 
61  // Construct search tree for all patch faces.
62  label sz = 0;
63  for (const label patchi : patchSet_)
64  {
65  const polyPatch& pp = mesh().boundaryMesh()[patchi];
66 
67  sz += pp.size();
68 
69  DebugInfo << " " << pp.name() << " size " << pp.size() << endl;
70  }
71 
72  labelList patchFaces(sz);
73  sz = 0;
74  treeBoundBox bb(boundBox::invertedBox);
75  for (const label patchi : patchSet_)
76  {
77  const polyPatch& pp = mesh().boundaryMesh()[patchi];
78 
79  forAll(pp, i)
80  {
81  patchFaces[sz++] = pp.start()+i;
82  }
83 
84  // Without reduction.
85  bb.add(pp.points(), pp.meshPoints());
86  }
87 
88  // Not very random
89  Random rndGen(123456);
90  // Make bb asymmetric just to avoid problems on symmetric meshes
91  bb = bb.extend(rndGen, 1e-4);
92 
93  // Make sure bb is 3D.
94  bb.min() -= point::uniform(ROOTVSMALL);
95  bb.max() += point::uniform(ROOTVSMALL);
96 
97 
98  indexedOctree<treeDataFace> patchTree
99  (
100  treeDataFace // all information needed to search faces
101  (
102  false, // do not cache bb
103  mesh(),
104  patchFaces // boundary faces only
105  ),
106  bb, // overall search domain
107  8, // maxLevel
108  10, // leafsize
109  3.0 // duplicity
110  );
111 
112  // Force calculation of face-diagonal decomposition
113  (void)mesh().tetBasePtIs();
114 
115 
116  // All the info for nearest. Construct to miss
117  List<mappedPatchBase::nearInfo> nearest(sampleCoords_.size());
118 
119  forAll(sampleCoords_, sampleI)
120  {
121  const point& sample = sampleCoords_[sampleI];
122 
123  pointIndexHit& nearInfo = nearest[sampleI].first();
124 
125  // Find the nearest locally
126  if (patchFaces.size())
127  {
128  nearInfo = patchTree.findNearest(sample, sqr(searchDist_));
129  }
130  else
131  {
132  nearInfo.setMiss();
133  }
134 
135 
136  // Fill in the distance field and the processor field
137  if (!nearInfo.hit())
138  {
139  nearest[sampleI].second().first() = Foam::sqr(GREAT);
140  nearest[sampleI].second().second() = Pstream::myProcNo();
141  }
142  else
143  {
144  // Set nearest to mesh face label
145  nearInfo.setIndex(patchFaces[nearInfo.index()]);
146 
147  nearest[sampleI].second().first() = magSqr
148  (
149  nearInfo.hitPoint()
150  - sample
151  );
152  nearest[sampleI].second().second() = Pstream::myProcNo();
153  }
154  }
155 
156 
157  // Find nearest.
158  Pstream::listCombineGather(nearest, mappedPatchBase::nearestEqOp());
160 
161 
162  if (debug && Pstream::master())
163  {
164  OFstream str
165  (
166  mesh().time().path()
167  / name()
168  + "_nearest.obj"
169  );
170  Info<< "Dumping mapping as lines from supplied points to"
171  << " nearest patch face to file " << str.name() << endl;
172 
173  label vertI = 0;
174 
175  forAll(nearest, i)
176  {
177  if (nearest[i].first().hit())
178  {
179  meshTools::writeOBJ(str, sampleCoords_[i]);
180  ++vertI;
181  meshTools::writeOBJ(str, nearest[i].first().hitPoint());
182  ++vertI;
183  str << "l " << vertI-1 << ' ' << vertI << nl;
184  }
185  }
186  }
187 
188 
189  // Store the sampling locations on the nearest processor
190  forAll(nearest, sampleI)
191  {
192  const pointIndexHit& nearInfo = nearest[sampleI].first();
193 
194  if (nearInfo.hit())
195  {
196  if (nearest[sampleI].second().second() == Pstream::myProcNo())
197  {
198  label facei = nearInfo.index();
199 
200  samplingPts.append(nearInfo.hitPoint());
201  samplingCells.append(mesh().faceOwner()[facei]);
202  samplingFaces.append(facei);
203  samplingSegments.append(0);
204  samplingCurveDist.append(1.0 * sampleI);
205  }
206  }
207  else
208  {
209  // No processor found point near enough. Mark with special value
210  // which is intercepted when interpolating
211  if (Pstream::master())
212  {
213  samplingPts.append(sampleCoords_[sampleI]);
214  samplingCells.append(-1);
215  samplingFaces.append(-1);
216  samplingSegments.append(0);
217  samplingCurveDist.append(1.0 * sampleI);
218  }
219  }
220  }
221 }
222 
223 
224 void Foam::patchCloudSet::genSamples()
225 {
226  // Storage for sample points
227  DynamicList<point> samplingPts;
228  DynamicList<label> samplingCells;
229  DynamicList<label> samplingFaces;
230  DynamicList<label> samplingSegments;
231  DynamicList<scalar> samplingCurveDist;
232 
233  calcSamples
234  (
235  samplingPts,
236  samplingCells,
237  samplingFaces,
238  samplingSegments,
239  samplingCurveDist
240  );
241 
242  samplingPts.shrink();
243  samplingCells.shrink();
244  samplingFaces.shrink();
245  samplingSegments.shrink();
246  samplingCurveDist.shrink();
247 
248  setSamples
249  (
250  samplingPts,
251  samplingCells,
252  samplingFaces,
253  samplingSegments,
254  samplingCurveDist
255  );
256 
257  if (debug)
258  {
259  write(Info);
260  }
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
265 
267 (
268  const word& name,
269  const polyMesh& mesh,
270  const meshSearch& searchEngine,
271  const word& axis,
272  const List<point>& sampleCoords,
273  const labelHashSet& patchSet,
274  const scalar searchDist
275 )
276 :
277  sampledSet(name, mesh, searchEngine, axis),
278  sampleCoords_(sampleCoords),
279  patchSet_(patchSet),
280  searchDist_(searchDist)
281 {
282  genSamples();
283 }
284 
285 
287 (
288  const word& name,
289  const polyMesh& mesh,
290  const meshSearch& searchEngine,
291  const dictionary& dict
292 )
293 :
294  sampledSet(name, mesh, searchEngine, dict),
295  sampleCoords_(dict.get<pointField>("points")),
296  patchSet_
297  (
298  mesh.boundaryMesh().patchSet(dict.get<wordRes>("patches"))
299  ),
300  searchDist_(dict.get<scalar>("maxDistance"))
301 {
302  genSamples();
303 }
304 
305 
306 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::sampledSet
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:83
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::coordSet::name
const word & name() const
Definition: coordSet.H:125
meshTools.H
patchCloudSet.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::meshSearch
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:60
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::uniform
static Vector< Cmpt > uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:164
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Foam::boundBox::invertedBox
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::patchCloudSet::patchCloudSet
patchCloudSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const word &axis, const List< point > &sampleCoords, const labelHashSet &patchSet, const scalar searchDist)
Construct from components.
Definition: patchCloudSet.C:267
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
treeDataFace.H
Foam::Field< vector >
treeBoundBox.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Time.H
Foam::sampledSet::mesh
const polyMesh & mesh() const
Definition: sampledSet.H:281
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::polyBoundaryMesh::patchSet
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
Definition: polyBoundaryMesh.C:864
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
rndGen
Random rndGen
Definition: createFields.H:23
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
sample
Minimal example by using system/controlDict.functions:
mappedPatchBase.H
Foam::polyMesh::tetBasePtIs
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:892