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-2022 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 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 template<class T>
41 struct isNotEqOp
42 {
43  void operator()(T& x, const T& y) const
44  {
45  const T unsetVal(-VGREAT*pTraits<T>::one);
46 
47  if (x != unsetVal)
48  {
49  // Keep x.
50 
51  // Note: should check for y != unsetVal but multiple sample cells
52  // already handled in read().
53  }
54  else
55  {
56  // x is not set. y might be.
57  x = y;
58  }
59  }
60 };
61 
62 } // End namespace Foam
63 
64 
65 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
66 
67 template<class GeoField>
69 Foam::probes::getOrLoadField(const word& fieldName) const
70 {
71  tmp<GeoField> tfield;
72 
73  if (loadFromFiles_)
74  {
75  tfield.reset
76  (
77  new GeoField
78  (
79  IOobject
80  (
81  fieldName,
82  mesh_.time().timeName(),
83  mesh_,
86  false
87  ),
88  mesh_
89  )
90  );
91  }
92  else
93  {
94  // Slightly paranoid here
95  tfield.cref(mesh_.cfindObject<GeoField>(fieldName));
96  }
97 
98  return tfield;
99 }
100 
101 
102 template<class Type>
104 (
105  const word& fieldName,
106  const Field<Type>& values
107 )
108 {
109  const MinMax<Type> limits(values);
110  const Type avgVal = average(values);
111 
112  this->setResult("average(" + fieldName + ")", avgVal);
113  this->setResult("min(" + fieldName + ")", limits.min());
114  this->setResult("max(" + fieldName + ")", limits.max());
115  this->setResult("size(" + fieldName + ")", values.size());
116 
117  if (verbose_)
118  {
119  Info<< name() << " : " << fieldName << nl
120  << " avg: " << avgVal << nl
121  << " min: " << limits.min() << nl
122  << " max: " << limits.max() << nl << nl;
123  }
124 }
125 
126 
127 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
128 
129 template<class Type>
130 void Foam::probes::writeValues
131 (
132  const word& fieldName,
133  const Field<Type>& values,
134  const scalar timeValue
135 )
136 {
137  if (Pstream::master())
138  {
139  const unsigned int width(IOstream::defaultPrecision() + 7);
140  OFstream& os = *probeFilePtrs_[fieldName];
141 
142  os << setw(width) << timeValue;
143 
144  forAll(values, probei)
145  {
146  if (includeOutOfBounds_ || processor_[probei] != -1)
147  {
148  os << ' ' << setw(width) << values[probei];
149  }
150  }
151  os << endl;
152  }
153 }
154 
155 
156 template<class GeoField>
157 void Foam::probes::performAction
158 (
159  const fieldGroup<GeoField>& fieldNames,
160  unsigned request
161 )
162 {
163  for (const word& fieldName : fieldNames)
164  {
165  tmp<GeoField> tfield = getOrLoadField<GeoField>(fieldName);
166 
167  if (tfield)
168  {
169  const auto& field = tfield();
170  const scalar timeValue = field.time().timeOutputValue();
171 
172  Field<typename GeoField::value_type> values(sample(field));
173 
174  this->storeResults(fieldName, values);
175  if (request & ACTION_WRITE)
176  {
177  this->writeValues(fieldName, values, timeValue);
178  }
179  }
180  }
181 }
182 
183 
184 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
185 
186 template<class Type>
188 Foam::probes::sample(const VolumeField<Type>& vField) const
189 {
190  const Type unsetVal(-VGREAT*pTraits<Type>::one);
191 
192  auto tvalues = tmp<Field<Type>>::New(Field<Type>(this->size(), unsetVal));
193  auto& values = tvalues.ref();
194 
195  if (fixedLocations_)
196  {
197  autoPtr<interpolation<Type>> interpPtr
198  (
199  interpolation<Type>::New(samplePointScheme_, vField)
200  );
201 
202  forAll(*this, probei)
203  {
204  if (elementList_[probei] >= 0)
205  {
206  const vector& position = operator[](probei);
207 
208  values[probei] = interpPtr().interpolate
209  (
210  position,
211  elementList_[probei],
212  -1
213  );
214  }
215  }
216  }
217  else
218  {
219  forAll(*this, probei)
220  {
221  if (elementList_[probei] >= 0)
222  {
223  values[probei] = vField[elementList_[probei]];
224  }
225  }
226  }
227 
228  Pstream::listCombineReduce(values, isNotEqOp<Type>());
230  return tvalues;
231 }
232 
233 
234 template<class Type>
236 Foam::probes::sample(const SurfaceField<Type>& sField) const
237 {
238  const Type unsetVal(-VGREAT*pTraits<Type>::one);
239 
240  auto tvalues = tmp<Field<Type>>::New(Field<Type>(this->size(), unsetVal));
241  auto& values = tvalues.ref();
242 
243  forAll(*this, probei)
244  {
245  if (faceList_[probei] >= 0)
246  {
247  values[probei] = sField[faceList_[probei]];
248  }
249  }
250 
251  Pstream::listCombineReduce(values, isNotEqOp<Type>());
253  return tvalues;
254 }
255 
256 
257 template<class Type>
259 Foam::probes::sample(const word& fieldName) const
260 {
261  return sample(mesh_.lookupObject<VolumeField<Type>>(fieldName));
262 }
263 
264 
265 template<class Type>
267 Foam::probes::sampleSurfaceField(const word& fieldName) const
268 {
269  return sample(mesh_.lookupObject<SurfaceField<Type>>(fieldName));
270 }
271 
272 
273 // ************************************************************************* //
Foam::surfaceFields.
rDeltaTY field()
void operator()(T &x, const T &y) const
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
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:487
A traits class, which is primarily used for primitives.
Definition: pTraits.H:50
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:196
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:416
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Ignore writing from objectRegistry::writeObject()
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:361
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.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
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:100
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
Generic templated field type.
Definition: Field.H:61
A class for handling words, derived from Foam::string.
Definition: word.H:63
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
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 time name of given scalar time formatted with the given precision.
Definition: Time.C:760
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)
Am I the master rank.
Definition: UPstream.H:672
const T & max() const noexcept
The max value (second)
Definition: MinMaxI.H:114
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
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
const fvMesh & mesh_
Reference to the fvMesh.
void reset(tmp< T > &&other) noexcept
Clear existing and transfer ownership.
Definition: tmpI.H:290
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.