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-2023 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 "dictionary.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "Time.H"
34 #include "IOmanip.H"
35 #include "mapPolyMesh.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(probes, 0);
43 
45  (
46  functionObject,
47  probes,
48  dictionary
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 void Foam::probes::createProbeFiles(const wordList& fieldNames)
56 {
57  // Open new output streams
58 
59  bool needsNewFiles = false;
60  for (const word& fieldName : fieldNames)
61  {
62  if (!probeFilePtrs_.found(fieldName))
63  {
64  needsNewFiles = true;
65  break;
66  }
67  }
68 
69  if (needsNewFiles && Pstream::master())
70  {
71  DebugInfo
72  << "Probing fields: " << fieldNames << nl
73  << "Probing locations: " << *this << nl
74  << endl;
75 
76  // Put in undecomposed case
77  // (Note: gives problems for distributed data running)
78 
79  fileName probeDir
80  (
83  / name()/mesh_.regionName()
84  // Use startTime as the instance for output files
86  );
87  probeDir.clean(); // Remove unneeded ".."
88 
89  // Create directory if needed
90  Foam::mkDir(probeDir);
91 
92  for (const word& fieldName : fieldNames)
93  {
94  if (probeFilePtrs_.found(fieldName))
95  {
96  // Safety
97  continue;
98  }
99 
100  auto osPtr = autoPtr<OFstream>::New(probeDir/fieldName);
101  auto& os = *osPtr;
102 
103  probeFilePtrs_.insert(fieldName, osPtr);
104 
105  DebugInfo<< "open probe stream: " << os.name() << endl;
106 
107  const unsigned int width(IOstream::defaultPrecision() + 7);
108  os << setf(ios_base::left);
109 
110  forAll(*this, probei)
111  {
112  os << "# Probe " << probei << ' ' << operator[](probei);
113 
114  if (processor_[probei] == -1)
115  {
116  os << " # Not Found";
117  }
118  // Only for patchProbes
119  else if (probei < patchIDList_.size())
120  {
121  const label patchi = patchIDList_[probei];
122  if (patchi != -1)
123  {
124  const polyBoundaryMesh& bm = mesh_.boundaryMesh();
125  if
126  (
127  patchi < bm.nNonProcessor()
128  || processor_[probei] == Pstream::myProcNo()
129  )
130  {
131  os << " at patch " << bm[patchi].name();
132  }
133  os << " with a distance of "
134  << mag(operator[](probei)-oldPoints_[probei])
135  << " m to the original point "
136  << oldPoints_[probei];
137  }
138  }
139 
140  os << nl;
141  }
142 
143  os << setw(width) << "# Time";
144 
145  forAll(*this, probei)
146  {
147  if (includeOutOfBounds_ || processor_[probei] != -1)
148  {
149  os << ' ' << setw(width) << probei;
150  }
151  }
152  os << endl;
153  }
154  }
155 }
156 
157 
158 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
159 
160 void Foam::probes::findElements(const fvMesh& mesh)
161 {
162  DebugInfo<< "probes: resetting sample locations" << endl;
163 
164  elementList_.resize_nocopy(pointField::size());
165  faceList_.resize_nocopy(pointField::size());
166  processor_.resize_nocopy(pointField::size());
167  processor_ = -1;
168 
169  forAll(*this, probei)
170  {
171  const point& location = (*this)[probei];
172 
173  const label celli = mesh.findCell(location);
174 
175  elementList_[probei] = celli;
176 
177  if (celli != -1)
178  {
179  const labelList& cellFaces = mesh.cells()[celli];
180  const vector& cellCentre = mesh.cellCentres()[celli];
181  scalar minDistance = GREAT;
182  label minFaceID = -1;
183  forAll(cellFaces, i)
184  {
185  label facei = cellFaces[i];
186  vector dist = mesh.faceCentres()[facei] - cellCentre;
187  if (mag(dist) < minDistance)
188  {
189  minDistance = mag(dist);
190  minFaceID = facei;
191  }
192  }
193  faceList_[probei] = minFaceID;
194  }
195  else
196  {
197  faceList_[probei] = -1;
198  }
199 
200  if (debug && (elementList_[probei] != -1 || faceList_[probei] != -1))
201  {
202  Pout<< "probes : found point " << location
203  << " in cell " << elementList_[probei]
204  << " and face " << faceList_[probei] << endl;
205  }
206  }
207 
208 
209  // Check if all probes have been found.
210  forAll(elementList_, probei)
211  {
212  const point& location = operator[](probei);
213  label celli = elementList_[probei];
214  label facei = faceList_[probei];
215 
216  processor_[probei] = (celli != -1 ? Pstream::myProcNo() : -1);
217 
218  // Check at least one processor with cell.
219  reduce(celli, maxOp<label>());
220  reduce(facei, maxOp<label>());
221  reduce(processor_[probei], maxOp<label>());
222 
223  if (celli == -1)
224  {
225  if (Pstream::master())
226  {
228  << "Did not find location " << location
229  << " in any cell. Skipping location." << endl;
230  }
231  }
232  else if (facei == -1)
233  {
234  if (Pstream::master())
235  {
237  << "Did not find location " << location
238  << " in any face. Skipping location." << endl;
239  }
240  }
241  else
242  {
243  // Make sure location not on two domains.
244  if (elementList_[probei] != -1 && elementList_[probei] != celli)
245  {
247  << "Location " << location
248  << " seems to be on multiple domains:"
249  << " cell " << elementList_[probei]
250  << " on my domain " << Pstream::myProcNo()
251  << " and cell " << celli << " on some other domain."
252  << nl
253  << "This might happen if the probe location is on"
254  << " a processor patch. Change the location slightly"
255  << " to prevent this." << endl;
256  }
257 
258  if (faceList_[probei] != -1 && faceList_[probei] != facei)
259  {
261  << "Location " << location
262  << " seems to be on multiple domains:"
263  << " cell " << faceList_[probei]
264  << " on my domain " << Pstream::myProcNo()
265  << " and face " << facei << " on some other domain."
266  << nl
267  << "This might happen if the probe location is on"
268  << " a processor patch. Change the location slightly"
269  << " to prevent this." << endl;
270  }
271  }
272  }
273 }
274 
275 
276 Foam::label Foam::probes::prepare(unsigned request)
277 {
278  // Prefilter on selection
279  HashTable<wordHashSet> selected =
280  (
281  loadFromFiles_
282  ? IOobjectList(mesh_, mesh_.time().timeName()).classes(fieldSelection_)
283  : mesh_.classes(fieldSelection_)
284  );
285 
286  // Classify and count fields
287  label nFields = 0;
288  do
289  {
290  #undef doLocalCode
291  #define doLocalCode(InputType, Target) \
292  { \
293  Target.clear(); /* Remove old values */ \
294  const auto iter = selected.cfind(InputType::typeName); \
295  if (iter.good()) \
296  { \
297  /* Add new (current) values */ \
298  Target.append(iter.val().sortedToc()); \
299  nFields += Target.size(); \
300  } \
301  }
302 
303  doLocalCode(volScalarField, scalarFields_);
304  doLocalCode(volVectorField, vectorFields_)
305  doLocalCode(volSphericalTensorField, sphericalTensorFields_);
306  doLocalCode(volSymmTensorField, symmTensorFields_);
307  doLocalCode(volTensorField, tensorFields_);
308 
309  doLocalCode(surfaceScalarField, surfaceScalarFields_);
310  doLocalCode(surfaceVectorField, surfaceVectorFields_);
311  doLocalCode(surfaceSphericalTensorField, surfaceSphericalTensorFields_);
312  doLocalCode(surfaceSymmTensorField, surfaceSymmTensorFields_);
313  doLocalCode(surfaceTensorField, surfaceTensorFields_);
314  #undef doLocalCode
315  }
316  while (false);
317 
318 
319  // Adjust file streams
320  if (Pstream::master())
321  {
322  wordHashSet currentFields(2*nFields);
323  currentFields.insert(scalarFields_);
324  currentFields.insert(vectorFields_);
325  currentFields.insert(sphericalTensorFields_);
326  currentFields.insert(symmTensorFields_);
327  currentFields.insert(tensorFields_);
328 
329  currentFields.insert(surfaceScalarFields_);
330  currentFields.insert(surfaceVectorFields_);
331  currentFields.insert(surfaceSphericalTensorFields_);
332  currentFields.insert(surfaceSymmTensorFields_);
333  currentFields.insert(surfaceTensorFields_);
334 
335  DebugInfo
336  << "Probing fields: " << currentFields << nl
337  << "Probing locations: " << *this << nl
338  << endl;
339 
340  // Close streams for fields that no longer exist
341  forAllIters(probeFilePtrs_, iter)
342  {
343  if (!currentFields.erase(iter.key()))
344  {
345  DebugInfo<< "close probe stream: " << iter()->name() << endl;
346 
347  probeFilePtrs_.remove(iter);
348  }
349  }
350 
351  if ((request & ACTION_WRITE) && !currentFields.empty())
352  {
353  createProbeFiles(currentFields.sortedToc());
354  }
355  }
356 
357  return nFields;
358 }
359 
360 
361 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
362 
363 Foam::probes::probes
364 (
365  const word& name,
366  const Time& runTime,
367  const dictionary& dict,
368  const bool loadFromFiles,
369  const bool readFields
370 )
371 :
372  functionObjects::fvMeshFunctionObject(name, runTime, dict),
373  pointField(),
374  loadFromFiles_(loadFromFiles),
375  fixedLocations_(true),
376  includeOutOfBounds_(true),
377  verbose_(false),
378  onExecute_(false),
379  fieldSelection_(),
380  samplePointScheme_("cell")
381 {
382  if (readFields)
383  {
385  }
386 }
387 
388 
389 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
390 
391 bool Foam::probes::verbose(const bool on) noexcept
392 {
393  bool old(verbose_);
394  verbose_ = on;
395  return old;
396 }
397 
398 
400 {
401  dict.readEntry("probeLocations", static_cast<pointField&>(*this));
402  dict.readEntry("fields", fieldSelection_);
403 
404  dict.readIfPresent("fixedLocations", fixedLocations_);
405  dict.readIfPresent("includeOutOfBounds", includeOutOfBounds_);
406 
407  verbose_ = dict.getOrDefault("verbose", false);
408  onExecute_ = dict.getOrDefault("sampleOnExecute", false);
409 
410  if (dict.readIfPresent("interpolationScheme", samplePointScheme_))
411  {
412  if (!fixedLocations_ && samplePointScheme_ != "cell")
413  {
415  << "Only cell interpolation can be applied when "
416  << "not using fixedLocations. InterpolationScheme "
417  << "entry will be ignored"
418  << endl;
419  }
420  }
421 
422  // Initialise cells to sample from supplied locations
423  findElements(mesh_);
424 
425  // Close old (ununsed) streams
426  prepare(ACTION_NONE);
427 
428  return true;
429 }
430 
431 
432 bool Foam::probes::performAction(unsigned request)
433 {
434  if (!pointField::empty() && request && prepare(request))
435  {
436  performAction(scalarFields_, request);
437  performAction(vectorFields_, request);
438  performAction(sphericalTensorFields_, request);
439  performAction(symmTensorFields_, request);
440  performAction(tensorFields_, request);
441 
442  performAction(surfaceScalarFields_, request);
443  performAction(surfaceVectorFields_, request);
444  performAction(surfaceSphericalTensorFields_, request);
445  performAction(surfaceSymmTensorFields_, request);
446  performAction(surfaceTensorFields_, request);
447  }
448  return true;
449 }
450 
451 
453 {
454  if (onExecute_)
455  {
456  return performAction(ACTION_ALL & ~ACTION_WRITE);
457  }
458 
459  return true;
460 }
461 
463 bool Foam::probes::write()
464 {
465  return performAction(ACTION_ALL);
466 }
467 
468 
469 void Foam::probes::updateMesh(const mapPolyMesh& mpm)
470 {
471  DebugInfo<< "probes: updateMesh" << endl;
472 
473  if (&mpm.mesh() != &mesh_)
474  {
475  return;
476  }
477 
478  if (fixedLocations_)
479  {
480  findElements(mesh_);
481  }
482  else
483  {
484  DebugInfo<< "probes: remapping sample locations" << endl;
485 
486  // 1. Update cells
487  {
488  DynamicList<label> elems(elementList_.size());
489 
490  const labelList& reverseMap = mpm.reverseCellMap();
491  forAll(elementList_, i)
492  {
493  label celli = elementList_[i];
494  if (celli != -1)
495  {
496  label newCelli = reverseMap[celli];
497  if (newCelli == -1)
498  {
499  // cell removed
500  }
501  else if (newCelli < -1)
502  {
503  // cell merged
504  elems.append(-newCelli - 2);
505  }
506  else
507  {
508  // valid new cell
509  elems.append(newCelli);
510  }
511  }
512  else
513  {
514  // Keep -1 elements so the size stays the same
515  elems.append(-1);
516  }
517  }
518 
519  elementList_.transfer(elems);
520  }
521 
522  // 2. Update faces
523  {
524  DynamicList<label> elems(faceList_.size());
525 
526  const labelList& reverseMap = mpm.reverseFaceMap();
527  for (const label facei : faceList_)
528  {
529  if (facei != -1)
530  {
531  label newFacei = reverseMap[facei];
532  if (newFacei == -1)
533  {
534  // face removed
535  }
536  else if (newFacei < -1)
537  {
538  // face merged
539  elems.append(-newFacei - 2);
540  }
541  else
542  {
543  // valid new face
544  elems.append(newFacei);
545  }
546  }
547  else
548  {
549  // Keep -1 elements
550  elems.append(-1);
551  }
552  }
554  faceList_.transfer(elems);
555  }
556  }
557 }
558 
559 
560 void Foam::probes::movePoints(const polyMesh& mesh)
561 {
562  DebugInfo<< "probes: movePoints" << endl;
563 
564  if (fixedLocations_ && &mesh == &mesh_)
565  {
566  findElements(mesh_);
567  }
568 }
569 
570 
571 // ************************************************************************* //
Foam::surfaceFields.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:85
const Type & value() const noexcept
Return const reference to value.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual bool write()
Sample and write.
Definition: probes.C:456
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
pointField oldPoints_
Original probes location (only used for patchProbes)
Definition: probes.H:286
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: probes.C:462
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
virtual bool execute()
Sample and store result if the sampleOnExecute is enabled.
Definition: probes.C:445
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:86
virtual const wordRes & fieldNames() const noexcept
Return names of fields to probe.
Definition: probes.H:402
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
#define doLocalCode(InputType, Target)
virtual dimensionedScalar startTime() const
Return start time.
Definition: Time.C:789
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
GeometricField< sphericalTensor, fvsPatchField, surfaceMesh > surfaceSphericalTensorField
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:423
labelList processor_
Processor holding the cell or face (-1 if point not found.
Definition: probes.H:271
const cellList & cells() const
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
Definition: volFieldsFwd.H:84
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
Macros for easy insertion into run-time selection tables.
const word & name() const noexcept
Return the name of this functionObject.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
label prepare(unsigned request)
Classify field types, close/open file streams.
Definition: probes.C:269
T & operator[](const label i)
Return element of UList.
Definition: UListI.H:361
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
dynamicFvMesh & mesh
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:614
GeometricField< symmTensor, fvsPatchField, surfaceMesh > surfaceSymmTensorField
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
fileName globalPath() const
Return global path for the case = rootPath/globalCaseName. Same as TimePaths::globalPath() ...
Definition: Time.H:512
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
GeometricField< tensor, fvsPatchField, surfaceMesh > surfaceTensorField
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:336
Smanip< ios_base::fmtflags > setf(const ios_base::fmtflags flags)
Definition: IOmanip.H:169
Vector< scalar > vector
Definition: vector.H:57
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:73
const vectorField & cellCentres() const
#define DebugInfo
Report an information message using Foam::Info.
const direction noexcept
Definition: Scalar.H:258
Istream and Ostream manipulators taking arguments.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
label size() const noexcept
The number of elements in the container.
Definition: UList.H:671
bool includeOutOfBounds_
Include probes that were not found (default: true)
Definition: probes.H:206
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:653
int debug
Static debugging option.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
HashPtrTable< OFstream > probeFilePtrs_
Current open files (non-empty on master only)
Definition: probes.H:276
virtual void movePoints(const polyMesh &)
Update for changes of mesh.
Definition: probes.C:553
const vectorField & faceCentres() const
List< word > wordList
List of word.
Definition: fileName.H:59
vector point
Point is a vector.
Definition: point.H:37
bool verbose(const bool on) noexcept
Enable/disable verbose output.
Definition: probes.C:384
#define WarningInFunction
Report a warning using Foam::Warning.
virtual void findElements(const fvMesh &mesh)
Find cells and faces containing probes.
Definition: probes.C:153
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:437
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
static word outputPrefix
Directory prefix.
virtual bool read(const dictionary &)
Read the probes.
Definition: probes.C:392
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1511
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
List< label > labelList
A List of labels.
Definition: List.H:62
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
labelList patchIDList_
Patch IDs on which the new probes are located (for patchProbes)
Definition: probes.H:281
const fvMesh & mesh_
Reference to the fvMesh.
static const word & regionName(const word &region)
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:811
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:617
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)