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-------------------------------------------------------------------------------
11License
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
41namespace Foam
42{
45}
46
47
48// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49
50void 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 - globally consistent
158 Pstream::listCombineAllGather(nearest, mappedPatchBase::nearestEqOp());
159
160
161 if (debug && Pstream::master())
162 {
163 OFstream str
164 (
165 mesh().time().path()
166 / name()
167 + "_nearest.obj"
168 );
169 Info<< "Dumping mapping as lines from supplied points to"
170 << " nearest patch face to file " << str.name() << endl;
171
172 label vertI = 0;
173
174 forAll(nearest, i)
175 {
176 if (nearest[i].first().hit())
177 {
178 meshTools::writeOBJ(str, sampleCoords_[i]);
179 ++vertI;
180 meshTools::writeOBJ(str, nearest[i].first().hitPoint());
181 ++vertI;
182 str << "l " << vertI-1 << ' ' << vertI << nl;
183 }
184 }
185 }
186
187
188 // Store the sampling locations on the nearest processor
189 forAll(nearest, sampleI)
190 {
191 const pointIndexHit& nearInfo = nearest[sampleI].first();
192
193 if (nearInfo.hit())
194 {
195 if (nearest[sampleI].second().second() == Pstream::myProcNo())
196 {
197 label facei = nearInfo.index();
198
199 samplingPts.append(nearInfo.hitPoint());
200 samplingCells.append(mesh().faceOwner()[facei]);
201 samplingFaces.append(facei);
202 samplingSegments.append(0);
203 samplingCurveDist.append(1.0 * sampleI);
204 }
205 }
206 else
207 {
208 // No processor found point near enough. Mark with special value
209 // which is intercepted when interpolating
210 if (Pstream::master())
211 {
212 samplingPts.append(sampleCoords_[sampleI]);
213 samplingCells.append(-1);
214 samplingFaces.append(-1);
215 samplingSegments.append(0);
216 samplingCurveDist.append(1.0 * sampleI);
217 }
218 }
219 }
220}
221
222
223void Foam::patchCloudSet::genSamples()
224{
225 // Storage for sample points
226 DynamicList<point> samplingPts;
227 DynamicList<label> samplingCells;
228 DynamicList<label> samplingFaces;
229 DynamicList<label> samplingSegments;
230 DynamicList<scalar> samplingCurveDist;
231
232 calcSamples
233 (
234 samplingPts,
235 samplingCells,
236 samplingFaces,
237 samplingSegments,
238 samplingCurveDist
239 );
240
241 samplingPts.shrink();
242 samplingCells.shrink();
243 samplingFaces.shrink();
244 samplingSegments.shrink();
245 samplingCurveDist.shrink();
246
247 setSamples
248 (
249 samplingPts,
250 samplingCells,
251 samplingFaces,
252 samplingSegments,
253 samplingCurveDist
254 );
255
256 if (debug)
257 {
258 write(Info);
259 }
260}
261
262
263// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
264
266(
267 const word& name,
268 const polyMesh& mesh,
269 const meshSearch& searchEngine,
270 const word& axis,
271 const List<point>& sampleCoords,
272 const labelHashSet& patchSet,
273 const scalar searchDist
274)
275:
276 sampledSet(name, mesh, searchEngine, axis),
277 sampleCoords_(sampleCoords),
278 patchSet_(patchSet),
279 searchDist_(searchDist)
280{
281 genSamples();
282}
283
284
286(
287 const word& name,
288 const polyMesh& mesh,
289 const meshSearch& searchEngine,
290 const dictionary& dict
291)
292:
293 sampledSet(name, mesh, searchEngine, dict),
294 sampleCoords_(dict.get<pointField>("points")),
295 patchSet_
296 (
297 mesh.boundaryMesh().patchSet(dict.get<wordRes>("patches"))
298 ),
299 searchDist_(dict.get<scalar>("maxDistance"))
300{
301 genSamples();
302}
303
304
305// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Minimal example by using system/controlDict.functions:
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
static void listCombineAllGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
T & first()
Return the first element of the list.
Definition: UListI.H:202
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
const word & name() const noexcept
The coord-set name.
Definition: coordSet.H:134
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:61
Like Foam::cloudSet but samples nearest patch face.
Definition: patchCloudSet.H:95
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:906
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
int myProcNo() const noexcept
Return processor number.
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:86
const polyMesh & mesh() const noexcept
Definition: sampledSet.H:319
splitCell * master() const
Definition: splitCell.H:113
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
dynamicFvMesh & mesh
#define DebugInfo
Report an information message using Foam::Info.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
dimensionedSymmTensor sqr(const dimensionedVector &dv)
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
messageStream Info
Information stream (stdout output on master, null elsewhere)
vector point
Point is a vector.
Definition: point.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
runTime write()
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Random rndGen
Definition: createFields.H:23