probes.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 "probes.H"
27 #include "volFields.H"
28 #include "dictionary.H"
29 #include "Time.H"
30 #include "IOmanip.H"
31 #include "mapPolyMesh.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(probes, 0);
39 
41  (
42  functionObject,
43  probes,
44  dictionary
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
52 {
53  if (debug)
54  {
55  Info<< "probes: resetting sample locations" << endl;
56  }
57 
58  elementList_.clear();
59  elementList_.setSize(size());
60 
61  faceList_.clear();
62  faceList_.setSize(size());
63 
64  forAll(*this, probei)
65  {
66  const vector& location = operator[](probei);
67 
68  const label celli = mesh.findCell(location);
69 
70  elementList_[probei] = celli;
71 
72  if (celli != -1)
73  {
74  const labelList& cellFaces = mesh.cells()[celli];
75  const vector& cellCentre = mesh.cellCentres()[celli];
76  scalar minDistance = GREAT;
77  label minFaceID = -1;
78  forAll(cellFaces, i)
79  {
80  label facei = cellFaces[i];
81  vector dist = mesh.faceCentres()[facei] - cellCentre;
82  if (mag(dist) < minDistance)
83  {
84  minDistance = mag(dist);
85  minFaceID = facei;
86  }
87  }
88  faceList_[probei] = minFaceID;
89  }
90  else
91  {
92  faceList_[probei] = -1;
93  }
94 
95  if (debug && (elementList_[probei] != -1 || faceList_[probei] != -1))
96  {
97  Pout<< "probes : found point " << location
98  << " in cell " << elementList_[probei]
99  << " and face " << faceList_[probei] << endl;
100  }
101  }
102 
103 
104  // Check if all probes have been found.
105  forAll(elementList_, probei)
106  {
107  const vector& location = operator[](probei);
108  label celli = elementList_[probei];
109  label facei = faceList_[probei];
110 
111  // Check at least one processor with cell.
112  reduce(celli, maxOp<label>());
113  reduce(facei, maxOp<label>());
114 
115  if (celli == -1)
116  {
117  if (Pstream::master())
118  {
120  << "Did not find location " << location
121  << " in any cell. Skipping location." << endl;
122  }
123  }
124  else if (facei == -1)
125  {
126  if (Pstream::master())
127  {
129  << "Did not find location " << location
130  << " in any face. Skipping location." << endl;
131  }
132  }
133  else
134  {
135  // Make sure location not on two domains.
136  if (elementList_[probei] != -1 && elementList_[probei] != celli)
137  {
139  << "Location " << location
140  << " seems to be on multiple domains:"
141  << " cell " << elementList_[probei]
142  << " on my domain " << Pstream::myProcNo()
143  << " and cell " << celli << " on some other domain."
144  << endl
145  << "This might happen if the probe location is on"
146  << " a processor patch. Change the location slightly"
147  << " to prevent this." << endl;
148  }
149 
150  if (faceList_[probei] != -1 && faceList_[probei] != facei)
151  {
153  << "Location " << location
154  << " seems to be on multiple domains:"
155  << " cell " << faceList_[probei]
156  << " on my domain " << Pstream::myProcNo()
157  << " and face " << facei << " on some other domain."
158  << endl
159  << "This might happen if the probe location is on"
160  << " a processor patch. Change the location slightly"
161  << " to prevent this." << endl;
162  }
163  }
164  }
165 }
166 
167 
169 {
170  const label nFields = classifyFields();
171 
172  // adjust file streams
173  if (Pstream::master())
174  {
175  wordHashSet currentFields;
176 
177  currentFields.insert(scalarFields_);
178  currentFields.insert(vectorFields_);
179  currentFields.insert(sphericalTensorFields_);
180  currentFields.insert(symmTensorFields_);
181  currentFields.insert(tensorFields_);
182 
183  currentFields.insert(surfaceScalarFields_);
184  currentFields.insert(surfaceVectorFields_);
185  currentFields.insert(surfaceSphericalTensorFields_);
186  currentFields.insert(surfaceSymmTensorFields_);
187  currentFields.insert(surfaceTensorFields_);
188 
189  if (debug)
190  {
191  Info<< "Probing fields: " << currentFields << nl
192  << "Probing locations: " << *this << nl
193  << endl;
194  }
195 
196 
197  fileName probeDir;
198  fileName probeSubDir = name();
199 
200  if (mesh_.name() != polyMesh::defaultRegion)
201  {
202  probeSubDir = probeSubDir/mesh_.name();
203  }
204  probeSubDir = "postProcessing"/probeSubDir/mesh_.time().timeName();
205 
206  if (Pstream::parRun())
207  {
208  // Put in undecomposed case
209  // (Note: gives problems for distributed data running)
210  probeDir = mesh_.time().path()/".."/probeSubDir;
211  }
212  else
213  {
214  probeDir = mesh_.time().path()/probeSubDir;
215  }
216 
217  // ignore known fields, close streams for fields that no longer exist
218  forAllIter(HashPtrTable<OFstream>, probeFilePtrs_, iter)
219  {
220  if (!currentFields.erase(iter.key()))
221  {
222  if (debug)
223  {
224  Info<< "close probe stream: " << iter()->name() << endl;
225  }
226 
227  delete probeFilePtrs_.remove(iter);
228  }
229  }
230 
231  // currentFields now just has the new fields - open streams for them
232  forAllConstIter(wordHashSet, currentFields, iter)
233  {
234  const word& fieldName = iter.key();
235 
236  // Create directory if does not exist.
237  mkDir(probeDir);
238 
239  OFstream* fPtr = new OFstream(probeDir/fieldName);
240 
241  OFstream& fout = *fPtr;
242 
243  if (debug)
244  {
245  Info<< "open probe stream: " << fout.name() << endl;
246  }
247 
248  probeFilePtrs_.insert(fieldName, fPtr);
249 
250  unsigned int w = IOstream::defaultPrecision() + 7;
251 
252  forAll(*this, probei)
253  {
254  fout<< "# Probe " << probei << ' ' << operator[](probei)
255  << endl;
256  }
257 
258  fout<< '#' << setw(IOstream::defaultPrecision() + 6)
259  << "Probe";
260 
261  forAll(*this, probei)
262  {
263  fout<< ' ' << setw(w) << probei;
264  }
265  fout<< endl;
266 
267  fout<< '#' << setw(IOstream::defaultPrecision() + 6)
268  << "Time" << endl;
269  }
270  }
271 
272  return nFields;
273 }
274 
275 
276 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
277 
278 Foam::probes::probes
279 (
280  const word& name,
281  const Time& t,
282  const dictionary& dict
283 )
284 :
285  functionObject(name),
286  pointField(0),
287  mesh_
288  (
289  refCast<const fvMesh>
290  (
292  (
294  )
295  )
296  ),
297  loadFromFiles_(false),
298  fieldSelection_(),
299  fixedLocations_(true),
300  interpolationScheme_("cell")
301 {
302  read(dict);
303 }
304 
305 
306 Foam::probes::probes
307 (
308  const word& name,
309  const objectRegistry& obr,
310  const dictionary& dict,
311  const bool loadFromFiles
312 )
313 :
314  functionObject(name),
315  pointField(0),
316  mesh_(refCast<const fvMesh>(obr)),
317  loadFromFiles_(loadFromFiles),
318  fieldSelection_(),
319  fixedLocations_(true),
320  interpolationScheme_("cell")
321 {
322  read(dict);
323 }
324 
325 
326 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
327 
329 {}
330 
331 
332 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
333 
335 {
336  dict.lookup("probeLocations") >> *this;
337  dict.lookup("fields") >> fieldSelection_;
338 
339  dict.readIfPresent("fixedLocations", fixedLocations_);
340  if (dict.readIfPresent("interpolationScheme", interpolationScheme_))
341  {
342  if (!fixedLocations_ && interpolationScheme_ != "cell")
343  {
345  << "Only cell interpolation can be applied when "
346  << "not using fixedLocations. InterpolationScheme "
347  << "entry will be ignored";
348  }
349  }
350 
351  // Initialise cells to sample from supplied locations
352  findElements(mesh_);
353 
354  prepare();
355 
356  return true;
357 }
358 
359 
361 {
362  return true;
363 }
364 
365 
367 {
368  if (size() && prepare())
369  {
370  sampleAndWrite(scalarFields_);
371  sampleAndWrite(vectorFields_);
372  sampleAndWrite(sphericalTensorFields_);
373  sampleAndWrite(symmTensorFields_);
374  sampleAndWrite(tensorFields_);
375 
376  sampleAndWriteSurfaceFields(surfaceScalarFields_);
377  sampleAndWriteSurfaceFields(surfaceVectorFields_);
378  sampleAndWriteSurfaceFields(surfaceSphericalTensorFields_);
379  sampleAndWriteSurfaceFields(surfaceSymmTensorFields_);
380  sampleAndWriteSurfaceFields(surfaceTensorFields_);
381  }
382 
383  return true;
384 }
385 
386 
388 {
389  DebugInfo<< "probes: updateMesh" << endl;
390 
391  if (&mpm.mesh() != &mesh_)
392  {
393  return;
394  }
395 
396  if (fixedLocations_)
397  {
398  findElements(mesh_);
399  }
400  else
401  {
402  if (debug)
403  {
404  Info<< "probes: remapping sample locations" << endl;
405  }
406 
407  // 1. Update cells
408  {
409  DynamicList<label> elems(elementList_.size());
410 
411  const labelList& reverseMap = mpm.reverseCellMap();
412  forAll(elementList_, i)
413  {
414  label celli = elementList_[i];
415  label newCelli = reverseMap[celli];
416  if (newCelli == -1)
417  {
418  // cell removed
419  }
420  else if (newCelli < -1)
421  {
422  // cell merged
423  elems.append(-newCelli - 2);
424  }
425  else
426  {
427  // valid new cell
428  elems.append(newCelli);
429  }
430  }
431 
432  elementList_.transfer(elems);
433  }
434 
435  // 2. Update faces
436  {
437  DynamicList<label> elems(faceList_.size());
438 
439  const labelList& reverseMap = mpm.reverseFaceMap();
440  forAll(faceList_, i)
441  {
442  label facei = faceList_[i];
443  label newFacei = reverseMap[facei];
444  if (newFacei == -1)
445  {
446  // face removed
447  }
448  else if (newFacei < -1)
449  {
450  // face merged
451  elems.append(-newFacei - 2);
452  }
453  else
454  {
455  // valid new face
456  elems.append(newFacei);
457  }
458  }
459 
460  faceList_.transfer(elems);
461  }
462  }
463 }
464 
465 
467 {
468  DebugInfo<< "probes: movePoints" << endl;
469 
470  if (fixedLocations_ && &mesh == &mesh_)
471  {
472  findElements(mesh_);
473  }
474 }
475 
476 
477 // ************************************************************************* //
A HashTable with keys but without contents.
Definition: HashSet.H:59
#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
A class for handling file names.
Definition: fileName.H:69
virtual bool write()
Sample and write.
Definition: probes.C:366
label prepare()
Classify field type and Open/close file streams,.
Definition: probes.C:168
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: probes.C:387
virtual bool execute()
Execute, currently does nothing.
Definition: probes.C:360
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
Output to file stream.
Definition: OFstream.H:81
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:417
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:411
A HashTable specialization for hashing pointers.
Definition: HashPtrTable.H:50
const vectorField & faceCentres() const
Abstract base-class for Time/database function objects.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
const cellList & cells() const
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Macros for easy insertion into run-time selection tables.
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Definition: HashTable.C:367
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
word name() const
Return file name (part beyond last /)
Definition: fileName.C:179
bool read(const char *, int32_t &)
Definition: int32IO.C:85
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dynamicFvMesh & mesh
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:495
A class for handling words, derived from string.
Definition: word.H:59
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
virtual void findElements(const fvMesh &)
Find cells and faces containing probes.
Definition: probes.C:51
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
const vectorField & cellCentres() const
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
#define DebugInfo
Report an information message using Foam::Info.
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
Istream and Ostream manipulators taking arguments.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:262
defineTypeNameAndDebug(combustionModel, 0)
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:295
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void movePoints(const polyMesh &)
Update for changes of mesh.
Definition: probes.C:466
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:393
#define WarningInFunction
Report a warning using Foam::Warning.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:526
virtual bool read(const dictionary &)
Read the probes.
Definition: probes.C:334
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1408
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Registry of regIOobjects.
virtual ~probes()
Destructor.
Definition: probes.C:328
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:249
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:118
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:357
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