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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "patchProbes.H"
27 #include "volFields.H"
28 #include "IOmanip.H"
29 // For 'nearInfo' helper class only
30 #include "mappedPatchBase.H"
31 #include "treeBoundBox.H"
32 #include "treeDataFace.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(patchProbes, 0);
39 }
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
44 {
45  (void)mesh.tetBasePtIs();
46 
47  const polyBoundaryMesh& bm = mesh.boundaryMesh();
48 
49  // All the info for nearest. Construct to miss
50  List<mappedPatchBase::nearInfo> nearest(this->size());
51 
52  const labelList patchIDs(bm.patchSet(patchNames_).sortedToc());
53 
54  label nFaces = 0;
55  forAll(patchIDs, i)
56  {
57  nFaces += bm[patchIDs[i]].size();
58  }
59 
60  if (nFaces > 0)
61  {
62  // Collect mesh faces and bounding box
63  labelList bndFaces(nFaces);
65 
66  nFaces = 0;
67  forAll(patchIDs, i)
68  {
69  const polyPatch& pp = bm[patchIDs[i]];
70  forAll(pp, i)
71  {
72  bndFaces[nFaces++] = pp.start()+i;
73  const face& f = pp[i];
74  forAll(f, fp)
75  {
76  const point& pt = pp.points()[f[fp]];
77  overallBb.min() = min(overallBb.min(), pt);
78  overallBb.max() = max(overallBb.max(), pt);
79  }
80  }
81  }
82 
83  Random rndGen(123456);
84  overallBb = overallBb.extend(rndGen, 1e-4);
85  overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
86  overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
87 
88  const indexedOctree<treeDataFace> boundaryTree
89  (
90  treeDataFace // all information needed to search faces
91  (
92  false, // do not cache bb
93  mesh,
94  bndFaces // patch faces only
95  ),
96  overallBb, // overall search domain
97  8, // maxLevel
98  10, // leafsize
99  3.0 // duplicity
100  );
101 
102 
103  forAll(probeLocations(), probeI)
104  {
105  const point sample = probeLocations()[probeI];
106 
107  scalar span = boundaryTree.bb().mag();
108 
109  pointIndexHit info = boundaryTree.findNearest
110  (
111  sample,
112  Foam::sqr(span)
113  );
114 
115  if (!info.hit())
116  {
117  info = boundaryTree.findNearest(sample, Foam::sqr(GREAT));
118  }
119 
120  label faceI = boundaryTree.shapes().faceLabels()[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.hitPoint();
142 
143  mappedPatchBase::nearInfo sampleInfo;
144 
145  sampleInfo.first() = pointIndexHit
146  (
147  true,
148  facePt,
149  faceI
150  );
151 
152  sampleInfo.second().first() = magSqr(facePt-sample);
153  sampleInfo.second().second() = Pstream::myProcNo();
154 
155  nearest[probeI]= sampleInfo;
156  }
157  }
158  }
159 
160 
161  // Find nearest.
164 
165 
166  // Update actual probe locations
167  forAll(nearest, sampleI)
168  {
169  operator[](sampleI) = nearest[sampleI].first().rawPoint();
170  }
171 
172 
173  if (debug)
174  {
175  Info<< "patchProbes::findElements" << " : " << endl;
176  forAll(nearest, sampleI)
177  {
178  label procI = nearest[sampleI].second().second();
179  label localI = nearest[sampleI].first().index();
180 
181  Info<< " " << sampleI << " coord:"<< operator[](sampleI)
182  << " found on processor:" << procI
183  << " in local face:" << localI
184  << " with location:" << nearest[sampleI].first().rawPoint()
185  << endl;
186  }
187  }
188 
189 
190  // Extract any local faces to sample
191  elementList_.setSize(nearest.size());
192  elementList_ = -1;
193  faceList_.setSize(nearest.size());
194  faceList_ = -1;
195 
196  forAll(nearest, sampleI)
197  {
198  if (nearest[sampleI].second().second() == Pstream::myProcNo())
199  {
200  // Store the face to sample
201  faceList_[sampleI] = nearest[sampleI].first().index();
202  }
203  }
204 }
205 
206 
208 {
209  if (!dict.readIfPresent("patches", patchNames_))
210  {
211  word patchName(dict.lookup("patchName"));
212  patchNames_ = wordReList(1, wordRe(patchName));
213  }
215 }
216 
217 
218 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
219 
221 (
222  const word& name,
223  const objectRegistry& obr,
224  const dictionary& dict,
225  const bool loadFromFiles,
226  const bool doFindElements
227 )
228 :
229  probes(name, obr, dict, loadFromFiles, false)
230 {
231  readDict(dict);
232 
233  if (doFindElements)
234  {
235  // Find the elements
236  findElements(mesh_);
237 
238  // Open the probe streams
239  prepare();
240  }
241 }
242 
243 
244 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
245 
247 {}
248 
249 
251 {
252  if (this->size() && prepare())
253  {
254  sampleAndWrite(scalarFields_);
255  sampleAndWrite(vectorFields_);
256  sampleAndWrite(sphericalTensorFields_);
257  sampleAndWrite(symmTensorFields_);
258  sampleAndWrite(tensorFields_);
259 
260  sampleAndWriteSurfaceFields(surfaceScalarFields_);
261  sampleAndWriteSurfaceFields(surfaceVectorFields_);
262  sampleAndWriteSurfaceFields(surfaceSphericalTensorFields_);
263  sampleAndWriteSurfaceFields(surfaceSymmTensorFields_);
264  sampleAndWriteSurfaceFields(surfaceTensorFields_);
265  }
266 }
267 
268 
270 {
271  readDict(dict);
272 
273  // Find the elements
274  findElements(mesh_);
275 
276  // Open the probe streams
277  prepare();
278 }
279 
280 
281 // ************************************************************************* //
Foam::patchProbes::findElements
virtual void findElements(const fvMesh &)
Find elements containing patchProbes.
Definition: patchProbes.C:43
volFields.H
Foam::patchProbes::patchNames_
wordReList patchNames_
Patches to sample.
Definition: patchProbes.H:101
Foam::Random
Simple random number generator.
Definition: Random.H:49
Foam::PointIndexHit::index
label index() const
Return index.
Definition: PointIndexHit.H:111
Foam::boundBox::mag
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::patchProbes::~patchProbes
virtual ~patchProbes()
Destructor.
Definition: patchProbes.C:246
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
Foam::Tuple2::second
const Type2 & second() const
Return second.
Definition: Tuple2.H:106
Foam::patchProbes::patchProbes
patchProbes(const patchProbes &)
Disallow default bitwise copy construct.
Foam::patchProbes::readDict
void readDict(const dictionary &dict)
Read dictionary settings.
Definition: patchProbes.C:207
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
Foam::Tuple2::first
const Type1 & first() const
Return first.
Definition: Tuple2.H:94
Foam::wordReList
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::PointIndexHit::hit
bool hit() const
Is there a hit.
Definition: PointIndexHit.H:105
Foam::indexedOctree::shapes
const Type & shapes() const
Reference to shape.
Definition: indexedOctree.H:449
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mappedPatchBase::nearestEqOp
Definition: mappedPatchBase.H:139
Foam::wordRe
A wordRe is a word, but can also have a regular expression for matching words.
Definition: wordRe.H:74
treeDataFace.H
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
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:59
Foam::treeBoundBox::extend
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
Definition: treeBoundBoxI.H:317
treeBoundBox.H
Foam::Info
messageStream Info
Foam::PointIndexHit::hitPoint
const Point & hitPoint() const
Return hit point.
Definition: PointIndexHit.H:117
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
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: treeDataTriSurface.H:47
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::probes::faceList_
labelList faceList_
Definition: probes.H:173
Foam::probes::elementList_
labelList elementList_
Definition: probes.H:170
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:703
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::probes
Set of locations to sample.
Definition: probes.H:102
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::HashTable::sortedToc
List< Key > sortedToc() const
Return the table of contents as a sorted list.
Definition: HashTable.C:216
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::indexedOctree::bb
const treeBoundBox & bb() const
Top bounding box.
Definition: indexedOctree.H:468
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:282
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::treeBoundBox::invertedBox
static const treeBoundBox invertedBox
As per boundBox::invertedBox, but with GREAT instead of VGREAT.
Definition: treeBoundBox.H:99
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
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:59
patchProbes.H
Foam::treeDataFace
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:58
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::patchProbes::read
virtual void read(const dictionary &)
Read.
Definition: patchProbes.C:269
Foam::indexedOctree::findNearest
void findNearest(const label nodeI, const linePointRef &ln, treeBoundBox &tightest, label &nearestShapeI, point &linePoint, point &nearestPoint, const FindNearestOp &fnOp) const
Find nearest point to line.
Definition: indexedOctree.C:590
Foam::patchProbes::write
virtual void write()
Public members.
Definition: patchProbes.C:250
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Tuple2
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:47
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:417
Foam::point
vector point
Point is a vector.
Definition: point.H:41
Foam::polyBoundaryMesh::patchSet
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
Definition: polyBoundaryMesh.C:750
Foam::probes::readDict
void readDict(const dictionary &dict)
Read dictionary settings.
Definition: probes.C:266
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
rndGen
cachedRandom rndGen(label(0), -1)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
mappedPatchBase.H