probesTemplates.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-2016 OpenFOAM Foundation
9  Copyright (C) 2017-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 "volFields.H"
31 #include "surfaceFields.H"
32 #include "IOmanip.H"
33 #include "interpolation.H"
34 #include "SpanStream.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 template<class T>
42 struct isNotEqOp
43 {
44  void operator()(T& x, const T& y) const
45  {
46  const T unsetVal(-VGREAT*pTraits<T>::one);
47 
48  if (x != unsetVal)
49  {
50  // Keep x.
51 
52  // Note: should check for y != unsetVal but multiple sample cells
53  // already handled in read().
54  }
55  else
56  {
57  // x is not set. y might be.
58  x = y;
59  }
60  }
61 };
62 
63 } // End namespace Foam
64 
65 
66 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
67 
68 template<class GeoField>
70 Foam::probes::getOrLoadField(const word& fieldName) const
71 {
72  tmp<GeoField> tfield;
73 
74  if (loadFromFiles_)
75  {
76  tfield.emplace
77  (
78  IOobject
79  (
80  fieldName,
81  mesh_.time().timeName(),
82  mesh_.thisDb(),
86  ),
87  mesh_
88  );
89  }
90  else
91  {
92  tfield.cref(mesh_.cfindObject<GeoField>(fieldName));
93  }
94 
95  return tfield;
96 }
97 
98 
99 template<class Type>
101 (
102  const word& fieldName,
103  const Field<Type>& values
104 )
105 {
106  const MinMax<Type> limits(values);
107  const Type avgVal = average(values);
108 
109  this->setResult("average(" + fieldName + ")", avgVal);
110  this->setResult("min(" + fieldName + ")", limits.min());
111  this->setResult("max(" + fieldName + ")", limits.max());
112  this->setResult("size(" + fieldName + ")", values.size());
113 
114  if (verbose_)
115  {
116  Info<< name() << " : " << fieldName << nl
117  << " avg: " << avgVal << nl
118  << " min: " << limits.min() << nl
119  << " max: " << limits.max() << nl << nl;
120  }
121 }
122 
123 
124 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
125 
126 template<class Type>
127 void Foam::probes::writeValues
128 (
129  const word& fieldName,
130  const Field<Type>& values,
131  const scalar timeValue
132 )
133 {
134  if (Pstream::master())
135  {
136  const unsigned int width(IOstream::defaultPrecision() + 7);
137  OFstream& os = *probeFilePtrs_[fieldName];
138 
139  os << setw(width) << timeValue;
140 
141  OCharStream buf;
142 
143  forAll(values, probei)
144  {
145  if (includeOutOfBounds_ || processor_[probei] != -1)
146  {
147  buf.rewind();
148  buf << values[probei];
149  os << ' ' << setw(width) << buf.str().data();
150  }
151  }
152  os << endl;
153  }
154 }
155 
156 
157 template<class GeoField>
158 void Foam::probes::performAction
159 (
160  const fieldGroup<GeoField>& fieldNames,
161  unsigned request
162 )
163 {
164  for (const word& fieldName : fieldNames)
165  {
166  tmp<GeoField> tfield = getOrLoadField<GeoField>(fieldName);
167 
168  if (tfield)
169  {
170  const auto& field = tfield();
171  const scalar timeValue = field.time().timeOutputValue();
172 
173  Field<typename GeoField::value_type> values(sample(field));
174 
175  this->storeResults(fieldName, values);
176  if (request & ACTION_WRITE)
177  {
178  this->writeValues(fieldName, values, timeValue);
179  }
180  }
181  }
182 }
183 
184 
185 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186 
187 template<class Type>
189 Foam::probes::sample(const VolumeField<Type>& vField) const
190 {
191  const Type unsetVal(-VGREAT*pTraits<Type>::one);
192 
193  auto tvalues = tmp<Field<Type>>::New(Field<Type>(this->size(), unsetVal));
194  auto& values = tvalues.ref();
195 
196  if (fixedLocations_)
197  {
198  autoPtr<interpolation<Type>> interpPtr
199  (
200  interpolation<Type>::New(samplePointScheme_, vField)
201  );
202 
203  forAll(*this, probei)
204  {
205  if (elementList_[probei] >= 0)
206  {
207  const vector& position = operator[](probei);
208 
209  values[probei] = interpPtr().interpolate
210  (
211  position,
212  elementList_[probei],
213  -1
214  );
215  }
216  }
217  }
218  else
219  {
220  forAll(*this, probei)
221  {
222  if (elementList_[probei] >= 0)
223  {
224  values[probei] = vField[elementList_[probei]];
225  }
226  }
227  }
228 
229  Pstream::listCombineReduce(values, isNotEqOp<Type>());
231  return tvalues;
232 }
233 
234 
235 template<class Type>
237 Foam::probes::sample(const SurfaceField<Type>& sField) const
238 {
239  const Type unsetVal(-VGREAT*pTraits<Type>::one);
240 
241  auto tvalues = tmp<Field<Type>>::New(Field<Type>(this->size(), unsetVal));
242  auto& values = tvalues.ref();
243 
244  forAll(*this, probei)
245  {
246  if (faceList_[probei] >= 0)
247  {
248  values[probei] = sField[faceList_[probei]];
249  }
250  }
251 
252  Pstream::listCombineReduce(values, isNotEqOp<Type>());
254  return tvalues;
255 }
256 
257 
258 template<class Type>
260 Foam::probes::sample(const word& fieldName) const
261 {
262  return sample(mesh_.lookupObject<VolumeField<Type>>(fieldName));
263 }
264 
265 
266 template<class Type>
268 Foam::probes::sampleSurfaceField(const word& fieldName) const
269 {
270  return sample(mesh_.lookupObject<SurfaceField<Type>>(fieldName));
271 }
272 
273 
274 // ************************************************************************* //
Foam::surfaceFields.
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1)
rDeltaTY field()
Input/output streams with (internal or external) character storage.
void operator()(T &x, const T &y) const
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
A min/max value pair with additional methods. In addition to conveniently storing values...
Definition: HashSet.H:72
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
tmp< GeoField > getOrLoadField(const word &fieldName) const
Get from registry or load from disk.
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition: tmpI.H:221
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:423
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
tmp< Field< Type > > sampleSurfaceField(const word &fieldName) const
Sample a surface field at all locations.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
scalar y
const T & min() const noexcept
The min value (first)
Definition: MinMaxI.H:107
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Generic templated field type.
Definition: Field.H:62
A class for handling words, derived from Foam::string.
Definition: word.H:63
Istream and Ostream manipulators taking arguments.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
OBJstream os(runTime.globalPath()/outputName)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
bool loadFromFiles_
Load fields from files (not from objectRegistry)
Definition: probes.H:193
void storeResults(const word &fieldName, const Field< Type > &values)
Store results: min/max/average/size.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
const T & max() const noexcept
The max value (second)
Definition: MinMaxI.H:121
tmp< Field< Type > > sample(const VolumeField< Type > &) const
Sample a volume field at all locations.
Abstract base class for volume field interpolation.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
A class for managing temporary objects.
Definition: HashPtrTable.H:50
T & emplace(Args &&... args)
Reset with emplace construction. Return reference to the new content.
Definition: tmpI.H:357
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
const fvMesh & mesh_
Reference to the fvMesh.
Do not request registration (bool: false)
Namespace for OpenFOAM.
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
After completion all processors have the same data.