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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-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 "probes.H"
30 #include "volFields.H"
31 #include "dictionary.H"
32 #include "Time.H"
33 #include "IOmanip.H"
34 #include "mapPolyMesh.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(probes, 0);
42 
44  (
45  functionObject,
46  probes,
47  dictionary
48  );
49 }
50 
51 
52 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
53 
55 {
56  DebugInfo<< "probes: resetting sample locations" << endl;
57 
59  elementList_.setSize(size());
60 
61  faceList_.clear();
62  faceList_.setSize(size());
63 
64  processor_.setSize(size());
65  processor_ = -1;
66 
67  forAll(*this, probei)
68  {
69  const vector& location = operator[](probei);
70 
71  const label celli = mesh.findCell(location);
72 
73  elementList_[probei] = celli;
74 
75  if (celli != -1)
76  {
77  const labelList& cellFaces = mesh.cells()[celli];
78  const vector& cellCentre = mesh.cellCentres()[celli];
79  scalar minDistance = GREAT;
80  label minFaceID = -1;
81  forAll(cellFaces, i)
82  {
83  label facei = cellFaces[i];
84  vector dist = mesh.faceCentres()[facei] - cellCentre;
85  if (mag(dist) < minDistance)
86  {
87  minDistance = mag(dist);
88  minFaceID = facei;
89  }
90  }
91  faceList_[probei] = minFaceID;
92  }
93  else
94  {
95  faceList_[probei] = -1;
96  }
97 
98  if (debug && (elementList_[probei] != -1 || faceList_[probei] != -1))
99  {
100  Pout<< "probes : found point " << location
101  << " in cell " << elementList_[probei]
102  << " and face " << faceList_[probei] << endl;
103  }
104  }
105 
106 
107  // Check if all probes have been found.
108  forAll(elementList_, probei)
109  {
110  const vector& location = operator[](probei);
111  label celli = elementList_[probei];
112  label facei = faceList_[probei];
113 
114  processor_[probei] = (celli != -1 ? Pstream::myProcNo() : -1);
115 
116  // Check at least one processor with cell.
117  reduce(celli, maxOp<label>());
118  reduce(facei, maxOp<label>());
119  reduce(processor_[probei], maxOp<label>());
120 
121  if (celli == -1)
122  {
123  if (Pstream::master())
124  {
126  << "Did not find location " << location
127  << " in any cell. Skipping location." << endl;
128  }
129  }
130  else if (facei == -1)
131  {
132  if (Pstream::master())
133  {
135  << "Did not find location " << location
136  << " in any face. Skipping location." << endl;
137  }
138  }
139  else
140  {
141  // Make sure location not on two domains.
142  if (elementList_[probei] != -1 && elementList_[probei] != celli)
143  {
145  << "Location " << location
146  << " seems to be on multiple domains:"
147  << " cell " << elementList_[probei]
148  << " on my domain " << Pstream::myProcNo()
149  << " and cell " << celli << " on some other domain."
150  << nl
151  << "This might happen if the probe location is on"
152  << " a processor patch. Change the location slightly"
153  << " to prevent this." << endl;
154  }
155 
156  if (faceList_[probei] != -1 && faceList_[probei] != facei)
157  {
159  << "Location " << location
160  << " seems to be on multiple domains:"
161  << " cell " << faceList_[probei]
162  << " on my domain " << Pstream::myProcNo()
163  << " and face " << facei << " on some other domain."
164  << nl
165  << "This might happen if the probe location is on"
166  << " a processor patch. Change the location slightly"
167  << " to prevent this." << endl;
168  }
169  }
170  }
171 }
172 
173 
175 {
176  const label nFields = classifyFields();
177 
178  // adjust file streams
179  if (Pstream::master())
180  {
181  wordHashSet currentFields;
182 
183  currentFields.insert(scalarFields_);
184  currentFields.insert(vectorFields_);
185  currentFields.insert(sphericalTensorFields_);
186  currentFields.insert(symmTensorFields_);
187  currentFields.insert(tensorFields_);
188 
189  currentFields.insert(surfaceScalarFields_);
190  currentFields.insert(surfaceVectorFields_);
191  currentFields.insert(surfaceSphericalTensorFields_);
192  currentFields.insert(surfaceSymmTensorFields_);
193  currentFields.insert(surfaceTensorFields_);
194 
195  DebugInfo
196  << "Probing fields: " << currentFields << nl
197  << "Probing locations: " << *this << nl
198  << endl;
199 
200 
201  fileName probeSubDir = name();
202 
203  if (mesh_.name() != polyMesh::defaultRegion)
204  {
205  probeSubDir = probeSubDir/mesh_.name();
206  }
207 
208  // Put in undecomposed case
209  // (Note: gives problems for distributed data running)
210 
211  fileName probeDir
212  (
213  mesh_.time().globalPath()
215  / probeSubDir
216  / mesh_.time().timeName()
217  );
218  probeDir.clean(); // Remove unneeded ".."
219 
220  // ignore known fields, close streams for fields that no longer exist
221  forAllIters(probeFilePtrs_, iter)
222  {
223  if (!currentFields.erase(iter.key()))
224  {
225  DebugInfo<< "close probe stream: " << iter()->name() << endl;
226 
227  probeFilePtrs_.remove(iter);
228  }
229  }
230 
231  // currentFields now just has the new fields - open streams for them
232  for (const word& fieldName : currentFields)
233  {
234  // Create directory if does not exist.
235  mkDir(probeDir);
236 
237  auto fPtr = autoPtr<OFstream>::New(probeDir/fieldName);
238  auto& fout = *fPtr;
239 
240  DebugInfo<< "open probe stream: " << fout.name() << endl;
241 
242  probeFilePtrs_.insert(fieldName, fPtr);
243 
244  unsigned int w = IOstream::defaultPrecision() + 7;
245 
246  forAll(*this, probei)
247  {
248  fout<< "# Probe " << probei << ' ' << operator[](probei);
249 
250  if (processor_[probei] == -1)
251  {
252  fout<< " # Not Found";
253  }
254  // Only for patchProbes
255  else if (probei < patchIDList_.size())
256  {
257  const label patchi = patchIDList_[probei];
258  if (patchi != -1)
259  {
260  const polyBoundaryMesh& bm = mesh_.boundaryMesh();
261  if
262  (
263  patchi < bm.nNonProcessor()
264  || processor_[probei] == Pstream::myProcNo()
265  )
266  {
267  fout<< " at patch " << bm[patchi].name();
268  }
269  fout<< " with a distance of "
270  << mag(operator[](probei)-oldPoints_[probei])
271  << " m to the original point "
272  << oldPoints_[probei];
273  }
274  }
275 
276  fout<< endl;
277  }
278 
279  fout<< '#' << setw(IOstream::defaultPrecision() + 6)
280  << "Probe";
281 
282  forAll(*this, probei)
283  {
284  if (includeOutOfBounds_ || processor_[probei] != -1)
285  {
286  fout<< ' ' << setw(w) << probei;
287  }
288  }
289  fout<< endl;
290 
291  fout<< '#' << setw(IOstream::defaultPrecision() + 6)
292  << "Time" << endl;
293  }
294  }
295 
296  return nFields;
297 }
298 
299 
300 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
301 
302 Foam::probes::probes
303 (
304  const word& name,
305  const Time& runTime,
306  const dictionary& dict,
307  const bool loadFromFiles,
308  const bool readFields
309 )
310 :
312  pointField(0),
313  mesh_
314  (
315  refCast<const fvMesh>
316  (
318  (
320  )
321  )
322  ),
323  loadFromFiles_(loadFromFiles),
324  fieldSelection_(),
325  fixedLocations_(true),
326  interpolationScheme_("cell"),
327  includeOutOfBounds_(true)
328 {
329  if (readFields)
330  {
331  read(dict);
332  }
333 }
334 
335 
336 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
337 
339 {
340  dict.readEntry("probeLocations", static_cast<pointField&>(*this));
341  dict.readEntry("fields", fieldSelection_);
342 
343  dict.readIfPresent("fixedLocations", fixedLocations_);
344  if (dict.readIfPresent("interpolationScheme", interpolationScheme_))
345  {
346  if (!fixedLocations_ && interpolationScheme_ != "cell")
347  {
349  << "Only cell interpolation can be applied when "
350  << "not using fixedLocations. InterpolationScheme "
351  << "entry will be ignored"
352  << endl;
353  }
354  }
355  dict.readIfPresent("includeOutOfBounds", includeOutOfBounds_);
356 
357  // Initialise cells to sample from supplied locations
358  findElements(mesh_);
359 
360  prepare();
361 
362  return true;
363 }
364 
365 
367 {
368  return true;
369 }
370 
371 
373 {
374  if (size() && prepare())
375  {
376  sampleAndWrite(scalarFields_);
377  sampleAndWrite(vectorFields_);
378  sampleAndWrite(sphericalTensorFields_);
379  sampleAndWrite(symmTensorFields_);
380  sampleAndWrite(tensorFields_);
381 
382  sampleAndWriteSurfaceFields(surfaceScalarFields_);
383  sampleAndWriteSurfaceFields(surfaceVectorFields_);
384  sampleAndWriteSurfaceFields(surfaceSphericalTensorFields_);
385  sampleAndWriteSurfaceFields(surfaceSymmTensorFields_);
386  sampleAndWriteSurfaceFields(surfaceTensorFields_);
387  }
388 
389  return true;
390 }
391 
392 
394 {
395  DebugInfo<< "probes: updateMesh" << endl;
396 
397  if (&mpm.mesh() != &mesh_)
398  {
399  return;
400  }
401 
402  if (fixedLocations_)
403  {
404  findElements(mesh_);
405  }
406  else
407  {
408  DebugInfo<< "probes: remapping sample locations" << endl;
409 
410  // 1. Update cells
411  {
412  DynamicList<label> elems(elementList_.size());
413 
414  const labelList& reverseMap = mpm.reverseCellMap();
415  forAll(elementList_, i)
416  {
417  label celli = elementList_[i];
418  if (celli != -1)
419  {
420  label newCelli = reverseMap[celli];
421  if (newCelli == -1)
422  {
423  // cell removed
424  }
425  else if (newCelli < -1)
426  {
427  // cell merged
428  elems.append(-newCelli - 2);
429  }
430  else
431  {
432  // valid new cell
433  elems.append(newCelli);
434  }
435  }
436  else
437  {
438  // Keep -1 elements so the size stays the same
439  elems.append(-1);
440  }
441  }
442 
443  elementList_.transfer(elems);
444  }
445 
446  // 2. Update faces
447  {
448  DynamicList<label> elems(faceList_.size());
449 
450  const labelList& reverseMap = mpm.reverseFaceMap();
451  for (const label facei : faceList_)
452  {
453  if (facei != -1)
454  {
455  label newFacei = reverseMap[facei];
456  if (newFacei == -1)
457  {
458  // face removed
459  }
460  else if (newFacei < -1)
461  {
462  // face merged
463  elems.append(-newFacei - 2);
464  }
465  else
466  {
467  // valid new face
468  elems.append(newFacei);
469  }
470  }
471  else
472  {
473  // Keep -1 elements
474  elems.append(-1);
475  }
476  }
477 
478  faceList_.transfer(elems);
479  }
480  }
481 }
482 
483 
485 {
486  DebugInfo<< "probes: movePoints" << endl;
487 
488  if (fixedLocations_ && &mesh == &mesh_)
489  {
490  findElements(mesh_);
491  }
492 }
493 
494 
495 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
volFields.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::maxOp
Definition: ops.H:223
Foam::probes::execute
virtual bool execute()
Execute, currently does nothing.
Definition: probes.C:366
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::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:318
Foam::DynamicList< label >
mapPolyMesh.H
Foam::fileName::name
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:199
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::functionObjects::stateFunctionObject
Base class for function objects, adding functionality to read/write state information (data required ...
Definition: stateFunctionObject.H:69
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::HashSet< word, Hash< word > >
probes.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
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::probes::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: probes.C:393
Foam::probes::read
virtual bool read(const dictionary &)
Read the probes.
Definition: probes.C:338
Foam::Field< vector >
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::probes::faceList_
labelList faceList_
Faces to be probed.
Definition: probes.H:181
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
Foam::probes::elementList_
labelList elementList_
Cells to be probed (obtained from the locations)
Definition: probes.H:178
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::functionObject::outputPrefix
static word outputPrefix
Directory prefix.
Definition: functionObject.H:376
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:223
Foam::probes::prepare
label prepare()
Classify field type and open/close file streams,.
Definition: probes.C:174
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
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::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::mapPolyMesh::reverseFaceMap
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:501
Foam::polyBoundaryMesh::nNonProcessor
label nNonProcessor() const
The number of patches before the first processor patch.
Definition: polyBoundaryMesh.C:535
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:342
Time.H
Foam::mapPolyMesh::reverseCellMap
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:532
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
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::Vector< scalar >
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dictionary.H
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
Foam::probes::write
virtual bool write()
Sample and write.
Definition: probes.C:372
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::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fileName::clean
static bool clean(std::string &str)
Definition: fileName.C:199
Foam::probes::findElements
virtual void findElements(const fvMesh &mesh)
Find cells and faces containing probes.
Definition: probes.C:54
Foam::probes::movePoints
virtual void movePoints(const polyMesh &)
Update for changes of mesh.
Definition: probes.C:484
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::mapPolyMesh::mesh
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:363
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
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::probes::processor_
labelList processor_
Processor holding the cell or face (-1 if point not found.
Definition: probes.H:185