patchProbes.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) 2016-2021 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 "patchProbes.H"
30#include "volFields.H"
31#include "IOmanip.H"
32#include "mappedPatchBase.H"
33#include "treeBoundBox.H"
34#include "treeDataFace.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
42
44 (
48 );
49}
50
51// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52
54{
55 (void)mesh.tetBasePtIs();
56
58
59 // All the info for nearest. Construct to miss
61
62 const labelList patchIDs(bm.patchSet(patchNames_).sortedToc());
63
64 label nFaces = 0;
65 forAll(patchIDs, i)
66 {
67 nFaces += bm[patchIDs[i]].size();
68 }
69
70 if (nFaces > 0)
71 {
72 // Collect mesh faces and bounding box
73 labelList bndFaces(nFaces);
75
76 nFaces = 0;
77 forAll(patchIDs, i)
78 {
79 const polyPatch& pp = bm[patchIDs[i]];
80 forAll(pp, i)
81 {
82 bndFaces[nFaces++] = pp.start()+i;
83 const face& f = pp[i];
84
85 // Without reduction.
86 overallBb.add(pp.points(), f);
87 }
88 }
89
90 Random rndGen(123456);
91 overallBb = overallBb.extend(rndGen, 1e-4);
92 overallBb.min() -= point::uniform(ROOTVSMALL);
93 overallBb.max() += point::uniform(ROOTVSMALL);
94
95
96 const indexedOctree<treeDataFace> boundaryTree
97 (
98 treeDataFace // all information needed to search faces
99 (
100 false, // do not cache bb
101 mesh,
102 bndFaces // patch faces only
103 ),
104 overallBb, // overall search domain
105 8, // maxLevel
106 10, // leafsize
107 3.0 // duplicity
108 );
109
110 forAll(probeLocations(), probei)
111 {
112 const point sample = probeLocations()[probei];
113
114 const scalar span = boundaryTree.bb().mag();
115
116 pointIndexHit info = boundaryTree.findNearest
117 (
118 sample,
119 Foam::sqr(span)
120 );
121
122 if (!info.hit())
123 {
124 info = boundaryTree.findNearest(sample, Foam::sqr(GREAT));
125 }
126
127 label facei = boundaryTree.shapes().faceLabels()[info.index()];
128
129 const label patchi = bm.whichPatch(facei);
130
131 if (isA<emptyPolyPatch>(bm[patchi]))
132 {
134 << " The sample point: " << sample
135 << " belongs to " << patchi
136 << " which is an empty patch. This is not permitted. "
137 << " This sample will not be included "
138 << endl;
139 }
140 else if (info.hit())
141 {
142 // Note: do we store the face centre or the actual nearest?
143 // We interpolate using the faceI only though (no
144 // interpolation) so it does not actually matter much, just for
145 // the location written to the header.
146
147 //const point& facePt = mesh.faceCentres()[faceI];
148 const point& facePt = info.hitPoint();
149
150 mappedPatchBase::nearInfo sampleInfo;
151
152 sampleInfo.first() = pointIndexHit
153 (
154 true,
155 facePt,
156 facei
157 );
158
159 sampleInfo.second().first() = magSqr(facePt - sample);
160 sampleInfo.second().second() = Pstream::myProcNo();
161
162 nearest[probei]= sampleInfo;
163 }
164 }
165 }
166
167
168 // Find nearest - globally consistent
170
171 oldPoints_.resize(this->size());
172
173 // Update actual probe locations and store old ones
174 forAll(nearest, samplei)
175 {
176 oldPoints_[samplei] = operator[](samplei);
177 operator[](samplei) = nearest[samplei].first().rawPoint();
178 }
179
180 if (debug)
181 {
183 forAll(nearest, samplei)
184 {
185 label proci = nearest[samplei].second().second();
186 label locali = nearest[samplei].first().index();
187
188 Info<< " " << samplei << " coord:"<< operator[](samplei)
189 << " found on processor:" << proci
190 << " in local face:" << locali
191 << " with location:" << nearest[samplei].first().rawPoint()
192 << endl;
193 }
194 }
195
196 // Extract any local faces to sample:
197 // - operator[] : actual point to sample (=nearest point on patch)
198 // - oldPoints_ : original provided point (might be anywhere in the mesh)
199 // - elementList_ : cells, not used
200 // - faceList_ : faces (now patch faces)
201 // - patchIDList_ : patch corresponding to faceList
202 // - processor_ : processor
204 elementList_ = -1;
205
206 faceList_.resize_nocopy(nearest.size());
207 faceList_ = -1;
208
209 processor_.resize_nocopy(nearest.size());
210 processor_ = -1;
211
213 patchIDList_ = -1;
214
215 forAll(nearest, sampleI)
216 {
217 processor_[sampleI] = nearest[sampleI].second().second();
218
219 if (nearest[sampleI].second().second() == Pstream::myProcNo())
220 {
221 // Store the face to sample
222 faceList_[sampleI] = nearest[sampleI].first().index();
223 const label facei = faceList_[sampleI];
224 if (facei != -1)
225 {
226 processor_[sampleI] = Pstream::myProcNo();
227 patchIDList_[sampleI] = bm.whichPatch(facei);
228 }
229 }
230 reduce(processor_[sampleI], maxOp<label>());
231 reduce(patchIDList_[sampleI], maxOp<label>());
232 }
233}
234
235
236// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
237
239(
240 const word& name,
241 const Time& runTime,
242 const dictionary& dict,
243 const bool loadFromFiles,
244 const bool readFields
245)
246:
247 probes(name, runTime, dict, loadFromFiles, false)
248{
249 if (readFields)
250 {
251 read(dict);
252 }
253}
254
255
256// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
257
258bool Foam::patchProbes::performAction(unsigned request)
259{
260 if (!pointField::empty() && request && prepare(request))
261 {
262 performAction(scalarFields_, request);
263 performAction(vectorFields_, request);
264 performAction(sphericalTensorFields_, request);
265 performAction(symmTensorFields_, request);
266 performAction(tensorFields_, request);
267
268 performAction(surfaceScalarFields_, request);
269 performAction(surfaceVectorFields_, request);
270 performAction(surfaceSphericalTensorFields_, request);
271 performAction(surfaceSymmTensorFields_, request);
272 performAction(surfaceTensorFields_, request);
273 }
274 return true;
275}
276
277
279{
280 if (onExecute_)
281 {
282 return performAction(ACTION_ALL & ~ACTION_WRITE);
283 }
284
285 return true;
286}
287
288
290{
291 return performAction(ACTION_ALL);
292}
293
294
296{
297 if (!dict.readIfPresent("patches", patchNames_))
298 {
299 patchNames_.resize(1);
300 patchNames_.first() = dict.get<word>("patch");
301 }
302
303 return probes::read(dict);
304}
305
306
307// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
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:
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
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
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:146
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:66
label index() const noexcept
Return the hit index.
bool hit() const noexcept
Is there a hit?
const point_type & hitPoint() const
Return hit point. Fatal if not hit.
const Field< point_type > & points() const noexcept
Return reference to global points.
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.
virtual bool read()
Re-read model coefficients if they have changed.
Random number generator.
Definition: Random.H:60
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: Tuple2.H:58
const T1 & first() const noexcept
Return first.
Definition: Tuple2.H:118
const T2 & second() const noexcept
Return second.
Definition: Tuple2.H:130
T & first()
Return the first element of the list.
Definition: UListI.H:202
label size() const noexcept
The number of elements in the UList.
Definition: UListI.H:420
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
T & operator[](const label i)
Return element of UList.
Definition: UListI.H:299
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:191
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:133
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Particle-size distribution model wherein random samples are drawn from the doubly-truncated uniform p...
Definition: uniform.H:164
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Abstract base-class for Time/database function objects.
static int debug
Flag to execute debug content.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Non-pointer based hierarchical recursive searching.
Definition: indexedOctree.H:74
const treeBoundBox & bb() const
Top bounding box.
const Type & shapes() const
Reference to shape.
Set of locations to sample at patches.
Definition: patchProbes.H:89
virtual void findElements(const fvMesh &mesh)
Find elements containing patchProbes.
Definition: patchProbes.C:53
wordRes patchNames_
Patches to sample.
Definition: patchProbes.H:95
virtual bool execute()
Sample and store result if the sampleOnExecute is enabled.
Definition: patchProbes.C:278
virtual bool write()
Sample and write.
Definition: patchProbes.C:289
virtual bool read(const dictionary &)
Read.
Definition: patchProbes.C:295
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
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.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:906
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
Set of locations to sample.
Definition: probes.H:166
virtual const pointField & probeLocations() const noexcept
Return locations to probe.
Definition: probes.H:341
labelList processor_
Processor holding the cell or face (-1 if point not found.
Definition: probes.H:242
pointField oldPoints_
Original probes location (only used for patchProbes)
Definition: probes.H:251
labelList elementList_
Cells to be probed (obtained from the locations)
Definition: probes.H:235
labelList patchIDList_
Patch IDs on which the new probes are located (for patchProbes)
Definition: probes.H:248
labelList faceList_
Faces to be probed.
Definition: probes.H:238
int myProcNo() const noexcept
Return processor number.
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
treeBoundBox extend(Random &rndGen, const scalar s) const
Return slightly wider bounding box.
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:60
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
dynamicFvMesh & mesh
engineTime & runTime
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
#define WarningInFunction
Report a warning using Foam::Warning.
#define InfoInFunction
Report an information message using Foam::Info.
Namespace for OpenFOAM.
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)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
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
labelList f(nPoints)
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