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 -------------------------------------------------------------------------------
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,
46  patchProbes,
47  dictionary
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 
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.
171 
172  oldPoints_.resize(this->size());
173 
174  // Update actual probe locations and store old ones
175  forAll(nearest, samplei)
176  {
177  oldPoints_[samplei] = operator[](samplei);
178  operator[](samplei) = nearest[samplei].first().rawPoint();
179  }
180 
181  if (debug)
182  {
183  InfoInFunction << nl;
184  forAll(nearest, samplei)
185  {
186  label proci = nearest[samplei].second().second();
187  label locali = nearest[samplei].first().index();
188 
189  Info<< " " << samplei << " coord:"<< operator[](samplei)
190  << " found on processor:" << proci
191  << " in local face:" << locali
192  << " with location:" << nearest[samplei].first().rawPoint()
193  << endl;
194  }
195  }
196 
197  // Extract any local faces to sample:
198  // - operator[] : actual point to sample (=nearest point on patch)
199  // - oldPoints_ : original provided point (might be anywhere in the mesh)
200  // - elementList_ : cells, not used
201  // - faceList_ : faces (now patch faces)
202  // - patchIDList_ : patch corresponding to faceList
203  // - processor_ : processor
204  elementList_.setSize(nearest.size());
205  elementList_ = -1;
206  faceList_.setSize(nearest.size());
207  faceList_ = -1;
208  processor_.setSize(nearest.size());
209  processor_ = -1;
210  patchIDList_.setSize(nearest.size());
211  patchIDList_ = -1;
212 
213  forAll(nearest, sampleI)
214  {
215  processor_[sampleI] = nearest[sampleI].second().second();
216 
217  if (nearest[sampleI].second().second() == Pstream::myProcNo())
218  {
219  // Store the face to sample
220  faceList_[sampleI] = nearest[sampleI].first().index();
221  const label facei = faceList_[sampleI];
222  if (facei != -1)
223  {
224  processor_[sampleI] = Pstream::myProcNo();
225  patchIDList_[sampleI] = bm.whichPatch(facei);
226  }
227  }
228  reduce(processor_[sampleI], maxOp<label>());
229  reduce(patchIDList_[sampleI], maxOp<label>());
230  }
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
235 
236 Foam::patchProbes::patchProbes
237 (
238  const word& name,
239  const Time& t,
240  const dictionary& dict,
241  const bool loadFromFiles,
242  const bool readFields
243 )
244 :
245  probes(name, t, dict, loadFromFiles, false)
246 {
247  if (readFields)
248  {
249  read(dict);
250  }
251 }
252 
253 
254 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
255 
257 {
258  if (this->size() && prepare())
259  {
260  sampleAndWrite(scalarFields_);
261  sampleAndWrite(vectorFields_);
262  sampleAndWrite(sphericalTensorFields_);
263  sampleAndWrite(symmTensorFields_);
264  sampleAndWrite(tensorFields_);
265 
266  sampleAndWriteSurfaceFields(surfaceScalarFields_);
267  sampleAndWriteSurfaceFields(surfaceVectorFields_);
268  sampleAndWriteSurfaceFields(surfaceSphericalTensorFields_);
269  sampleAndWriteSurfaceFields(surfaceSymmTensorFields_);
270  sampleAndWriteSurfaceFields(surfaceTensorFields_);
271  }
272 
273  return true;
274 }
275 
276 
278 {
279  if (!dict.readIfPresent("patches", patchNames_))
280  {
281  patchNames_.resize(1);
282  patchNames_.first() = dict.get<word>("patch");
283  }
284 
285  return probes::read(dict);
286 }
287 
288 
289 // ************************************************************************* //
Foam::PrimitivePatch::points
const Field< point_type > & points() const noexcept
Return reference to global points.
Definition: PrimitivePatch.H:299
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::patchProbes::findElements
virtual void findElements(const fvMesh &)
Find elements containing patchProbes.
Definition: patchProbes.C:53
volFields.H
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::maxOp
Definition: ops.H:223
Foam::probes::patchIDList_
labelList patchIDList_
Patch IDs on which the new probes are located.
Definition: probes.H:193
Foam::treeBoundBox::extend
treeBoundBox extend(Random &rndGen, const scalar s) const
Return slightly wider bounding box.
Definition: treeBoundBoxI.H:325
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
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::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::patchProbes::read
virtual bool read(const dictionary &)
Read.
Definition: patchProbes.C:277
Foam::boundBox::invertedBox
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:86
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::indexedOctree::shapes
const Type & shapes() const
Reference to shape.
Definition: indexedOctree.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::PointIndexHit::hitPoint
const point_type & hitPoint() const
Return hit point. Fatal if not hit.
Definition: PointIndexHit.H:154
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::mappedPatchBase::nearestEqOp
Definition: mappedPatchBase.H:149
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
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::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:52
Foam::PointIndexHit::hit
bool hit() const noexcept
Is there a hit?
Definition: PointIndexHit.H:130
Foam::probes::read
virtual bool read(const dictionary &)
Read the probes.
Definition: probes.C:338
treeBoundBox.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:50
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::probes::faceList_
labelList faceList_
Faces to be probed.
Definition: probes.H:181
Foam::probes::elementList_
labelList elementList_
Cells to be probed (obtained from the locations)
Definition: probes.H:178
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:812
Foam::probes
Set of locations to sample.
Definition: probes.H:110
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
Foam::probes::oldPoints_
pointField oldPoints_
Original probes location (only used for patchProbes)
Definition: probes.H:196
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:361
Foam::indexedOctree::bb
const treeBoundBox & bb() const
Top bounding box.
Definition: indexedOctree.H:463
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
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
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::functionObject::debug
static int debug
Flag to execute debug content.
Definition: functionObject.H:367
Foam::probes::probeLocations
virtual const pointField & probeLocations() const
Return locations to probe.
Definition: probes.H:279
f
labelList f(nPoints)
Foam::Vector< scalar >
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::PointIndexHit::index
label index() const noexcept
Return the hit index.
Definition: PointIndexHit.H:136
patchProbes.H
Foam::treeDataFace
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:59
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
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
Foam::functionObjects::readFields
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:155
Foam::Tuple2::second
const T2 & second() const noexcept
Return second.
Definition: Tuple2.H:130
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::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::Tuple2
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:60
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::patchProbes::patchNames_
wordRes patchNames_
Patches to sample.
Definition: patchProbes.H:101
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::patchProbes::write
virtual bool write()
Public members.
Definition: patchProbes.C:256
Foam::Tuple2::first
const T1 & first() const noexcept
Return first.
Definition: Tuple2.H:118
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405
Foam::boundBox::add
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:191
sample
Minimal example by using system/controlDict.functions:
mappedPatchBase.H
Foam::probes::processor_
labelList processor_
Processor holding the cell or face (-1 if point not found.
Definition: probes.H:185