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-2016 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 #include "mappedPatchBase.H"
30 #include "treeBoundBox.H"
31 #include "treeDataFace.H"
33 
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(patchProbes, 0);
40 
42  (
43  functionObject,
44  patchProbes,
45  dictionary
46  );
47 }
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
52 {
53  (void)mesh.tetBasePtIs();
54 
55  const polyBoundaryMesh& bm = mesh.boundaryMesh();
56 
57  label patchi = bm.findPatchID(patchName_);
58 
59  if (patchi == -1)
60  {
62  << " Unknown patch name "
63  << patchName_ << endl
64  << exit(FatalError);
65  }
66 
67  // All the info for nearest. Construct to miss
68  List<mappedPatchBase::nearInfo> nearest(this->size());
69 
70  const polyPatch& pp = bm[patchi];
71 
72  if (pp.size() > 0)
73  {
74  labelList bndFaces(pp.size());
75  forAll(bndFaces, i)
76  {
77  bndFaces[i] = pp.start() + i;
78  }
79 
80  treeBoundBox overallBb(pp.points());
81  Random rndGen(123456);
82  overallBb = overallBb.extend(rndGen, 1e-4);
83  overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
84  overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
85 
86  const indexedOctree<treeDataFace> boundaryTree
87  (
88  treeDataFace // all information needed to search faces
89  (
90  false, // do not cache bb
91  mesh,
92  bndFaces // patch faces only
93  ),
94  overallBb, // overall search domain
95  8, // maxLevel
96  10, // leafsize
97  3.0 // duplicity
98  );
99 
100 
101  forAll(probeLocations(), probei)
102  {
103  const point sample = probeLocations()[probei];
104 
105  scalar span = boundaryTree.bb().mag();
106 
107  pointIndexHit info = boundaryTree.findNearest
108  (
109  sample,
110  Foam::sqr(span)
111  );
112 
113  if (!info.hit())
114  {
115  info = boundaryTree.findNearest
116  (
117  sample,
118  Foam::sqr(GREAT)
119  );
120  }
121 
122  label facei = boundaryTree.shapes().faceLabels()[info.index()];
123 
124  const label patchi = bm.whichPatch(facei);
125 
126  if (isA<emptyPolyPatch>(bm[patchi]))
127  {
129  << " The sample point: " << sample
130  << " belongs to " << patchi
131  << " which is an empty patch. This is not permitted. "
132  << " This sample will not be included "
133  << endl;
134  }
135  else
136  {
137  const point& fc = mesh.faceCentres()[facei];
138 
139  mappedPatchBase::nearInfo sampleInfo;
140 
141  sampleInfo.first() = pointIndexHit
142  (
143  true,
144  fc,
145  facei
146  );
147 
148  sampleInfo.second().first() = magSqr(fc-sample);
149  sampleInfo.second().second() = Pstream::myProcNo();
150 
151  nearest[probei]= sampleInfo;
152  }
153  }
154  }
155 
156 
157  // Find nearest.
160 
161  if (debug)
162  {
163  InfoInFunction << endl;
164  forAll(nearest, sampleI)
165  {
166  label proci = nearest[sampleI].second().second();
167  label localI = nearest[sampleI].first().index();
168 
169  Info<< " " << sampleI << " coord:"<< operator[](sampleI)
170  << " found on processor:" << proci
171  << " in local cell/face:" << localI
172  << " with fc:" << nearest[sampleI].first().rawPoint() << endl;
173  }
174  }
175 
176 
177  // Extract any local faces to sample
178  elementList_.setSize(nearest.size(), -1);
179 
180  forAll(nearest, sampleI)
181  {
182  if (nearest[sampleI].second().second() == Pstream::myProcNo())
183  {
184  // Store the face to sample
185  elementList_[sampleI] = nearest[sampleI].first().index();
186  }
187  }
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
192 
193 Foam::patchProbes::patchProbes
194 (
195  const word& name,
196  const Time& t,
197  const dictionary& dict
198 )
199 :
200  probes(name, t, dict)
201 {
202  // When constructing probes above it will have called the
203  // probes::findElements (since the virtual mechanism not yet operating).
204  // Not easy to workaround (apart from feeding through flag into constructor)
205  // so clear out any cells found for now.
206  elementList_.clear();
207  faceList_.clear();
208 
209  read(dict);
210 }
211 
212 
213 Foam::patchProbes::patchProbes
214 (
215  const word& name,
216  const objectRegistry& obr,
217  const dictionary& dict,
218  const bool loadFromFiles
219 )
220 :
221  probes(name, obr, dict)
222 {
223  // When constructing probes above it will have called the
224  // probes::findElements (since the virtual mechanism not yet operating).
225  // Not easy to workaround (apart from feeding through flag into constructor)
226  // so clear out any cells found for now.
227  elementList_.clear();
228  faceList_.clear();
229 
230  read(dict);
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
235 
237 {}
238 
239 
241 {
242  if (this->size() && prepare())
243  {
244  sampleAndWrite(scalarFields_);
245  sampleAndWrite(vectorFields_);
246  sampleAndWrite(sphericalTensorFields_);
247  sampleAndWrite(symmTensorFields_);
248  sampleAndWrite(tensorFields_);
249 
250  sampleAndWriteSurfaceFields(surfaceScalarFields_);
251  sampleAndWriteSurfaceFields(surfaceVectorFields_);
252  sampleAndWriteSurfaceFields(surfaceSphericalTensorFields_);
253  sampleAndWriteSurfaceFields(surfaceSymmTensorFields_);
254  sampleAndWriteSurfaceFields(surfaceTensorFields_);
255  }
256 
257  return true;
258 }
259 
260 
262 {
263  dict.lookup("patchName") >> patchName_;
264  return probes::read(dict);
265 }
266 
267 
268 // ************************************************************************* //
cachedRandom rndGen(label(0),-1)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const Type1 & first() const
Return first.
Definition: Tuple2.H:94
const double e
Elementary charge.
Definition: doubleFloat.H:78
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:47
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:58
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
bool hit() const
Is there a hit.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual ~patchProbes()
Destructor.
Definition: patchProbes.C:236
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:417
Set of locations to sample.
Definition: probes.H:68
virtual bool read(const dictionary &)
Read.
Definition: patchProbes.C:261
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const Field< PointType > & points() const
Return reference to global points.
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
const vectorField & faceCentres() const
T & first()
Return the first element of the list.
Definition: UListI.H:114
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Macros for easy insertion into run-time selection tables.
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dynamicFvMesh & mesh
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:356
A class for handling words, derived from string.
Definition: word.H:59
Simple random number generator.
Definition: Random.H:49
Foam::polyBoundaryMesh.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Istream and Ostream manipulators taking arguments.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const treeBoundBox & bb() const
Top bounding box.
defineTypeNameAndDebug(combustionModel, 0)
const labelList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:818
const Type2 & second() const
Return second.
Definition: Tuple2.H:106
virtual void findElements(const fvMesh &)
Find elements containing patchProbes.
Definition: patchProbes.C:51
label patchi
vector point
Point is a vector.
Definition: point.H:41
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:47
#define WarningInFunction
Report a warning using Foam::Warning.
const Type & shapes() const
Reference to shape.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label index() const
Return index.
virtual bool write()
Public members.
Definition: patchProbes.C:240
virtual bool read(const dictionary &)
Read the probes.
Definition: probes.C:334
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
messageStream Info
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
label findPatchID(const word &patchName) const
Find patch index given a name.
const labelList & faceLabels() const
Definition: treeDataFace.H:179
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
#define InfoInFunction
Report an information message using Foam::Info.