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