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-2018 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 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(patchProbes, 0);
43 
45  (
46  functionObject,
47  patchProbes,
48  dictionary
49  );
50 }
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
55 {
56  (void)mesh.tetBasePtIs();
57 
58  const polyBoundaryMesh& bm = mesh.boundaryMesh();
59 
60  // All the info for nearest. Construct to miss
61  List<mappedPatchBase::nearInfo> nearest(this->size());
62 
63  const labelList patchIDs(bm.patchSet(patchNames_).sortedToc());
64 
65  label nFaces = 0;
66  forAll(patchIDs, i)
67  {
68  nFaces += bm[patchIDs[i]].size();
69  }
70 
71  if (nFaces > 0)
72  {
73  // Collect mesh faces and bounding box
74  labelList bndFaces(nFaces);
76 
77  nFaces = 0;
78  forAll(patchIDs, i)
79  {
80  const polyPatch& pp = bm[patchIDs[i]];
81  forAll(pp, i)
82  {
83  bndFaces[nFaces++] = pp.start()+i;
84  const face& f = pp[i];
85 
86  // Without reduction.
87  overallBb.add(pp.points(), f);
88  }
89  }
90 
91  Random rndGen(123456);
92  overallBb = overallBb.extend(rndGen, 1e-4);
93  overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
94  overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
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 
111  forAll(probeLocations(), probei)
112  {
113  const point sample = probeLocations()[probei];
114 
115  scalar span = boundaryTree.bb().mag();
116 
117  pointIndexHit info = boundaryTree.findNearest
118  (
119  sample,
120  Foam::sqr(span)
121  );
122 
123  if (!info.hit())
124  {
125  info = boundaryTree.findNearest(sample, Foam::sqr(GREAT));
126  }
127 
128  label facei = boundaryTree.shapes().faceLabels()[info.index()];
129 
130  const label patchi = bm.whichPatch(facei);
131 
132  if (isA<emptyPolyPatch>(bm[patchi]))
133  {
135  << " The sample point: " << sample
136  << " belongs to " << patchi
137  << " which is an empty patch. This is not permitted. "
138  << " This sample will not be included "
139  << endl;
140  }
141  else if (info.hit())
142  {
143  // Note: do we store the face centre or the actual nearest?
144  // We interpolate using the faceI only though (no
145  // interpolation) so it does not actually matter much, just for
146  // the location written to the header.
147 
148  //const point& facePt = mesh.faceCentres()[faceI];
149  const point& facePt = info.hitPoint();
150 
151  mappedPatchBase::nearInfo sampleInfo;
152 
153  sampleInfo.first() = pointIndexHit
154  (
155  true,
156  facePt,
157  facei
158  );
159 
160  sampleInfo.second().first() = magSqr(facePt - sample);
161  sampleInfo.second().second() = Pstream::myProcNo();
162 
163  nearest[probei]= sampleInfo;
164  }
165  }
166  }
167 
168 
169  // Find nearest.
172 
173 
174  // Update actual probe locations
175  forAll(nearest, samplei)
176  {
177  operator[](samplei) = nearest[samplei].first().rawPoint();
178  }
179 
180 
181  if (debug)
182  {
183  InfoInFunction << endl;
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 
198  // Extract any local faces to sample
199  elementList_.setSize(nearest.size());
200  elementList_ = -1;
201  faceList_.setSize(nearest.size());
202  faceList_ = -1;
203  processor_.setSize(nearest.size());
204  processor_ = -1;
205 
206  processor_.setSize(size());
207  processor_ = -1;
208 
209  forAll(nearest, sampleI)
210  {
211  processor_[sampleI] = nearest[sampleI].second().second();
212  if (nearest[sampleI].second().second() == Pstream::myProcNo())
213  {
214  // Store the face to sample
215  faceList_[sampleI] = nearest[sampleI].first().index();
216  label facei = faceList_[sampleI];
217  processor_[sampleI] = (facei != -1 ? Pstream::myProcNo() : -1);
218  }
219  reduce(processor_[sampleI], maxOp<label>());
220  }
221 }
222 
223 
224 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
225 
226 Foam::patchProbes::patchProbes
227 (
228  const word& name,
229  const Time& t,
230  const dictionary& dict,
231  const bool loadFromFiles,
232  const bool readFields
233 )
234 :
235  probes(name, t, dict, loadFromFiles, false)
236 {
237  if (readFields)
238  {
239  read(dict);
240  }
241 }
242 
243 
244 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
245 
247 {
248  if (this->size() && prepare())
249  {
250  sampleAndWrite(scalarFields_);
251  sampleAndWrite(vectorFields_);
252  sampleAndWrite(sphericalTensorFields_);
253  sampleAndWrite(symmTensorFields_);
254  sampleAndWrite(tensorFields_);
255 
256  sampleAndWriteSurfaceFields(surfaceScalarFields_);
257  sampleAndWriteSurfaceFields(surfaceVectorFields_);
258  sampleAndWriteSurfaceFields(surfaceSphericalTensorFields_);
259  sampleAndWriteSurfaceFields(surfaceSymmTensorFields_);
260  sampleAndWriteSurfaceFields(surfaceTensorFields_);
261  }
262 
263  return true;
264 }
265 
266 
268 {
269  if (!dict.readIfPresent("patches", patchNames_))
270  {
271  patchNames_.resize(1);
272  patchNames_.first() = dict.get<word>("patch");
273  }
274 
275  return probes::read(dict);
276 }
277 
278 
279 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::patchProbes::findElements
virtual void findElements(const fvMesh &)
Find elements containing patchProbes.
Definition: patchProbes.C:54
volFields.H
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::maxOp
Definition: ops.H:223
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:316
Foam::PointIndexHit::index
label index() const
Return index.
Definition: PointIndexHit.H:113
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatch.H:300
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::patchProbes::read
virtual bool read(const dictionary &)
Read.
Definition: patchProbes.C:267
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:87
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::PointIndexHit::hit
bool hit() const
Is there a hit.
Definition: PointIndexHit.H:107
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:337
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:81
Foam::mappedPatchBase::nearestEqOp
Definition: mappedPatchBase.H:143
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:290
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:55
Foam::probes::read
virtual bool read(const dictionary &)
Read the probes.
Definition: probes.C:317
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
treeBoundBox.H
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::PointIndexHit::hitPoint
const Point & hitPoint() const
Return hit point.
Definition: PointIndexHit.H:119
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::patchProbes::sample
tmp< Field< Type > > sample(const GeometricField< Type, fvPatchField, volMesh > &) const
Sample a volume field at all locations.
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:805
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:121
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:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::HashTable::sortedToc
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:136
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:311
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:444
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
Definition: pointIndexHit.H:45
Foam::functionObject::debug
static int debug
Definition: functionObject.H:253
Foam::probes::probeLocations
virtual const pointField & probeLocations() const
Return locations to probe.
Definition: probes.H:271
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: HashTable.H:102
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:857
Foam::functionObjects::readFields
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:108
Foam::Tuple2::second
const T2 & second() const noexcept
Return second.
Definition: Tuple2.H:130
rndGen
Random rndGen
Definition: createFields.H:23
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
Foam::Tuple2
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: Tuple2.H:57
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::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
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:246
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:294
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:417
Foam::boundBox::add
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:191
mappedPatchBase.H
Foam::probes::processor_
labelList processor_
Processor holding the cell or face (-1 if point not found.
Definition: probes.H:185