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 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 class isNotEqOp
42 {
43 public:
44 
45  void operator()(T& x, const T& y) const
46  {
47  const T unsetVal(-VGREAT*pTraits<T>::one);
48 
49  if (x != unsetVal)
50  {
51  // Keep x.
52 
53  // Note: should check for y != unsetVal but multiple sample cells
54  // already handled in read().
55  }
56  else
57  {
58  // x is not set. y might be.
59  x = y;
60  }
61  }
62 };
63 
64 }
65 
66 
67 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
68 
69 template<class Type>
70 void Foam::probes::sampleAndWrite
71 (
72  const GeometricField<Type, fvPatchField, volMesh>& vField
73 )
74 {
75  Field<Type> values(sample(vField));
76 
77  if (Pstream::master())
78  {
79  unsigned int w = IOstream::defaultPrecision() + 7;
80  OFstream& os = *probeFilePtrs_[vField.name()];
81 
82  os << setw(w) << vField.time().timeOutputValue();
83 
84  forAll(values, probei)
85  {
86  if (includeOutOfBounds_ || processor_[probei] != -1)
87  {
88  os << ' ' << setw(w) << values[probei];
89  }
90  }
91  os << endl;
92  }
93 }
94 
95 
96 template<class Type>
97 void Foam::probes::sampleAndWrite
98 (
99  const GeometricField<Type, fvsPatchField, surfaceMesh>& sField
100 )
101 {
102  Field<Type> values(sample(sField));
103 
104  if (Pstream::master())
105  {
106  unsigned int w = IOstream::defaultPrecision() + 7;
107  OFstream& os = *probeFilePtrs_[sField.name()];
108 
109  os << setw(w) << sField.time().timeOutputValue();
110 
111  forAll(values, probei)
112  {
113  if (includeOutOfBounds_ || processor_[probei] != -1)
114  {
115  os << ' ' << setw(w) << values[probei];
116  }
117  }
118  os << endl;
119  }
120 }
121 
122 
123 template<class Type>
124 void Foam::probes::sampleAndWrite(const fieldGroup<Type>& fields)
125 {
126  forAll(fields, fieldi)
127  {
128  if (loadFromFiles_)
129  {
130  sampleAndWrite
131  (
132  GeometricField<Type, fvPatchField, volMesh>
133  (
134  IOobject
135  (
136  fields[fieldi],
137  mesh_.time().timeName(),
138  mesh_,
141  false
142  ),
143  mesh_
144  )
145  );
146  }
147  else
148  {
149  objectRegistry::const_iterator iter = mesh_.find(fields[fieldi]);
150 
151  if
152  (
153  iter.found()
154  && iter()->type()
155  == GeometricField<Type, fvPatchField, volMesh>::typeName
156  )
157  {
158  sampleAndWrite
159  (
160  mesh_.lookupObject
161  <GeometricField<Type, fvPatchField, volMesh>>
162  (
163  fields[fieldi]
164  )
165  );
166  }
167  }
168  }
169 }
170 
171 
172 template<class Type>
173 void Foam::probes::sampleAndWriteSurfaceFields(const fieldGroup<Type>& fields)
174 {
175  forAll(fields, fieldi)
176  {
177  if (loadFromFiles_)
178  {
179  sampleAndWrite
180  (
181  GeometricField<Type, fvsPatchField, surfaceMesh>
182  (
183  IOobject
184  (
185  fields[fieldi],
186  mesh_.time().timeName(),
187  mesh_,
190  false
191  ),
192  mesh_
193  )
194  );
195  }
196  else
197  {
198  objectRegistry::const_iterator iter = mesh_.find(fields[fieldi]);
199 
200  if
201  (
202  iter.found()
203  && iter()->type()
204  == GeometricField<Type, fvsPatchField, surfaceMesh>::typeName
205  )
206  {
207  sampleAndWrite
208  (
209  mesh_.lookupObject
210  <GeometricField<Type, fvsPatchField, surfaceMesh>>
211  (
212  fields[fieldi]
213  )
214  );
215  }
216  }
217  }
218 }
219 
220 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221 
222 template<class Type>
225 (
227 ) const
228 {
229  const Type unsetVal(-VGREAT*pTraits<Type>::one);
230 
231  tmp<Field<Type>> tValues
232  (
233  new Field<Type>(this->size(), unsetVal)
234  );
235 
236  Field<Type>& values = tValues.ref();
237 
238  if (fixedLocations_)
239  {
240  autoPtr<interpolation<Type>> interpolator
241  (
242  interpolation<Type>::New(interpolationScheme_, vField)
243  );
244 
245  forAll(*this, probei)
246  {
247  if (elementList_[probei] >= 0)
248  {
249  const vector& position = operator[](probei);
250 
251  values[probei] = interpolator().interpolate
252  (
253  position,
254  elementList_[probei],
255  -1
256  );
257  }
258  }
259  }
260  else
261  {
262  forAll(*this, probei)
263  {
264  if (elementList_[probei] >= 0)
265  {
266  values[probei] = vField[elementList_[probei]];
267  }
268  }
269  }
270 
273 
274  return tValues;
275 }
276 
277 
278 template<class Type>
280 Foam::probes::sample(const word& fieldName) const
281 {
282  return sample
283  (
285  (
286  fieldName
287  )
288  );
289 }
290 
291 
292 template<class Type>
295 (
297 ) const
298 {
299  const Type unsetVal(-VGREAT*pTraits<Type>::one);
300 
301  tmp<Field<Type>> tValues
302  (
303  new Field<Type>(this->size(), unsetVal)
304  );
305 
306  Field<Type>& values = tValues.ref();
307 
308  forAll(*this, probei)
309  {
310  if (faceList_[probei] >= 0)
311  {
312  values[probei] = sField[faceList_[probei]];
313  }
314  }
315 
318 
319  return tValues;
320 }
321 
322 
323 template<class Type>
326 {
327  return sample
328  (
330  (
331  fieldName
332  )
333  );
334 }
335 
336 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
volFields.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Foam::isNotEqOp
Definition: probesTemplates.C:41
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
surfaceFields.H
Foam::surfaceFields.
interpolation.H
Foam::isNotEqOp::operator()
void operator()(T &x, const T &y) const
Definition: probesTemplates.C:45
probes.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::probes::sample
tmp< Field< Type > > sample(const GeometricField< Type, fvPatchField, volMesh > &) const
Sample a volume field at all locations.
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::probes::probeFilePtrs_
HashPtrTable< OFstream > probeFilePtrs_
Current open files.
Definition: probes.H:188
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::probes::includeOutOfBounds_
bool includeOutOfBounds_
Include probes that were not found.
Definition: probes.H:158
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::interpolation
Abstract base class for interpolation.
Definition: mappedPatchFieldBase.H:96
Foam::probes::sampleSurfaceFields
tmp< Field< Type > > sampleSurfaceFields(const word &fieldName) const
Sample a single scalar field on all sample locations.
os
OBJstream os(runTime.globalPath()/outputName)
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:342
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::Vector< scalar >
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::GeometricField< Type, fvPatchField, volMesh >
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
y
scalar y
Definition: LISASMDCalcMethod1.H:14
sample
Minimal example by using system/controlDict.functions:
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::probes::processor_
labelList processor_
Processor holding the cell or face (-1 if point not found.
Definition: probes.H:185