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-2022 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 "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 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(patchProbes, 0);
42 
44  (
45  functionObject,
48  );
49 }
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
54 {
55  (void)mesh.tetBasePtIs();
56 
57  const polyBoundaryMesh& bm = mesh.boundaryMesh();
58 
59  // All the info for nearest. Construct to miss
60  List<mappedPatchBase::nearInfo> nearest(this->size());
61 
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);
74  treeBoundBox overallBb;
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.inflate(rndGen, 1e-4, ROOTVSMALL);
92 
93 
94  const indexedOctree<treeDataFace> boundaryTree
95  (
96  treeDataFace(mesh, bndFaces), // patch faces only
97 
98  overallBb, // overall search domain
99  8, // maxLevel
100  10, // leafsize
101  3.0 // duplicity
102  );
103 
104  forAll(probeLocations(), probei)
105  {
106  const auto& treeData = boundaryTree.shapes();
107  const point sample = probeLocations()[probei];
108 
109  pointIndexHit info = boundaryTree.findNearest
110  (
111  sample,
112  Foam::sqr(boundaryTree.bb().mag())
113  );
114 
115  if (!info.hit())
116  {
117  info = boundaryTree.findNearest(sample, Foam::sqr(GREAT));
118  }
119 
120  const label facei = treeData.objectIndex(info.index());
121 
122  const label patchi = bm.whichPatch(facei);
123 
124  if (isA<emptyPolyPatch>(bm[patchi]))
125  {
127  << " The sample point: " << sample
128  << " belongs to " << patchi
129  << " which is an empty patch. This is not permitted. "
130  << " This sample will not be included "
131  << endl;
132  }
133  else if (info.hit())
134  {
135  // Note: do we store the face centre or the actual nearest?
136  // We interpolate using the faceI only though (no
137  // interpolation) so it does not actually matter much, just for
138  // the location written to the header.
139 
140  //const point& facePt = mesh.faceCentres()[faceI];
141  const point& facePt = info.point();
142 
143  mappedPatchBase::nearInfo sampleInfo;
144 
145  sampleInfo.first() = pointIndexHit(true, facePt, facei);
146 
147  sampleInfo.second().first() = facePt.distSqr(sample);
148  sampleInfo.second().second() = Pstream::myProcNo();
149 
150  nearest[probei] = sampleInfo;
151  }
152  }
153  }
154 
155 
156  // Find nearest - globally consistent
157  Pstream::listCombineReduce(nearest, mappedPatchBase::nearestEqOp());
158 
159  oldPoints_.resize(this->size());
160 
161  // Update actual probe locations and store old ones
162  forAll(nearest, samplei)
163  {
164  oldPoints_[samplei] = operator[](samplei);
165  operator[](samplei) = nearest[samplei].first().point();
166  }
167 
168  if (debug)
169  {
170  InfoInFunction << nl;
171  forAll(nearest, samplei)
172  {
173  label proci = nearest[samplei].second().second();
174  label locali = nearest[samplei].first().index();
175 
176  Info<< " " << samplei << " coord:"<< operator[](samplei)
177  << " found on processor:" << proci
178  << " in local face:" << locali
179  << " with location:" << nearest[samplei].first().point()
180  << endl;
181  }
182  }
183 
184  // Extract any local faces to sample:
185  // - operator[] : actual point to sample (=nearest point on patch)
186  // - oldPoints_ : original provided point (might be anywhere in the mesh)
187  // - elementList_ : cells, not used
188  // - faceList_ : faces (now patch faces)
189  // - patchIDList_ : patch corresponding to faceList
190  // - processor_ : processor
191  elementList_.resize_nocopy(nearest.size());
192  elementList_ = -1;
193 
194  faceList_.resize_nocopy(nearest.size());
195  faceList_ = -1;
196 
197  processor_.resize_nocopy(nearest.size());
198  processor_ = -1;
199 
200  patchIDList_.resize_nocopy(nearest.size());
201  patchIDList_ = -1;
202 
203  forAll(nearest, sampleI)
204  {
205  processor_[sampleI] = nearest[sampleI].second().second();
206 
207  if (nearest[sampleI].second().second() == Pstream::myProcNo())
208  {
209  // Store the face to sample
210  faceList_[sampleI] = nearest[sampleI].first().index();
211  const label facei = faceList_[sampleI];
212  if (facei != -1)
213  {
214  processor_[sampleI] = Pstream::myProcNo();
215  patchIDList_[sampleI] = bm.whichPatch(facei);
216  }
217  }
218  reduce(processor_[sampleI], maxOp<label>());
219  reduce(patchIDList_[sampleI], maxOp<label>());
220  }
221 }
222 
223 
224 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
225 
226 Foam::patchProbes::patchProbes
227 (
228  const word& name,
229  const Time& runTime,
230  const dictionary& dict,
231  const bool loadFromFiles,
232  const bool readFields
233 )
234 :
235  probes(name, runTime, dict, loadFromFiles, false)
236 {
237  if (readFields)
238  {
239  read(dict);
240  }
241 }
242 
243 
244 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
245 
246 bool Foam::patchProbes::performAction(unsigned request)
247 {
248  if (!pointField::empty() && request && prepare(request))
249  {
250  performAction(scalarFields_, request);
251  performAction(vectorFields_, request);
252  performAction(sphericalTensorFields_, request);
253  performAction(symmTensorFields_, request);
254  performAction(tensorFields_, request);
255 
256  performAction(surfaceScalarFields_, request);
257  performAction(surfaceVectorFields_, request);
258  performAction(surfaceSphericalTensorFields_, request);
259  performAction(surfaceSymmTensorFields_, request);
260  performAction(surfaceTensorFields_, request);
261  }
262  return true;
263 }
264 
265 
267 {
268  if (onExecute_)
269  {
270  return performAction(ACTION_ALL & ~ACTION_WRITE);
271  }
272 
273  return true;
274 }
275 
278 {
279  return performAction(ACTION_ALL);
280 }
281 
282 
284 {
285  if (!dict.readIfPresent("patches", patchNames_))
286  {
287  patchNames_.resize(1);
288  patchNames_.first() = dict.get<word>("patch");
289  }
290 
291  return probes::read(dict);
292 }
293 
294 
295 // ************************************************************************* //
const labelList patchIDs(pbm.indices(polyPatchNames, true))
dictionary dict
virtual void findElements(const fvMesh &mesh)
Find elements containing patchProbes.
Definition: patchProbes.C:46
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
pointField oldPoints_
Original probes location (only used for patchProbes)
Definition: probes.H:286
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
labelHashSet patchSet(const UList< wordRe > &select, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
Random rndGen
Definition: createFields.H:23
virtual bool read(const dictionary &)
Read.
Definition: patchProbes.C:276
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
wordRes patchNames_
Patches to sample.
Definition: patchProbes.H:93
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
labelList processor_
Processor holding the cell or face (-1 if point not found.
Definition: probes.H:271
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:175
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
T & operator[](const label i)
Return element of UList.
Definition: UListI.H:361
labelList faceList_
Faces to be probed.
Definition: probes.H:264
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
virtual const pointField & probeLocations() const noexcept
Return locations to probe.
Definition: probes.H:410
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
virtual bool execute()
Sample and store result if the sampleOnExecute is enabled.
Definition: patchProbes.C:259
A class for handling words, derived from Foam::string.
Definition: word.H:63
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
Set of locations to sample at patches.
Definition: patchProbes.H:82
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
Istream and Ostream manipulators taking arguments.
label size() const noexcept
The number of elements in the container.
Definition: UList.H:671
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
tmp< Field< Type > > sample(const VolumeField< Type > &) const
Sample a volume field at all locations.
Tuple2< pointIndexHit, Tuple2< scalar, label > > nearInfo
Helper class for finding nearest.
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:139
labelList elementList_
Cells to be probed (obtained from the locations)
Definition: probes.H:259
virtual bool write()
Sample and write.
Definition: patchProbes.C:270
virtual bool read(const dictionary &)
Read the probes.
Definition: probes.C:392
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< label > labelList
A List of labels.
Definition: List.H:62
labelList patchIDList_
Patch IDs on which the new probes are located (for patchProbes)
Definition: probes.H:281
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static int debug
Flag to execute debug content.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
After completion all processors have the same data.
#define InfoInFunction
Report an information message using Foam::Info.