probesTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
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 "surfaceFields.H"
29 #include "IOmanip.H"
30 #include "interpolation.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 template<class T>
38 class isNotEqOp
39 {
40 public:
41 
42  void operator()(T& x, const T& y) const
43  {
44  const T unsetVal(-vGreat*pTraits<T>::one);
45 
46  if (x != unsetVal)
47  {
48  // Keep x.
49 
50  // Note: should check for y != unsetVal but multiple sample cells
51  // already handled in read().
52  }
53  else
54  {
55  // x is not set. y might be.
56  x = y;
57  }
58  }
59 };
60 
61 }
62 
63 
64 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
65 
66 template<class Type>
67 void Foam::probes::sampleAndWrite
68 (
70 )
71 {
72  Field<Type> values(sample(vField));
73 
74  if (Pstream::master())
75  {
76  unsigned int w = IOstream::defaultPrecision() + 7;
77  OFstream& os = *probeFilePtrs_[vField.name()];
78 
79  os << setw(w) << vField.time().timeToUserTime(vField.time().value());
80 
81  forAll(values, probei)
82  {
83  os << ' ' << setw(w) << values[probei];
84  }
85  os << endl;
86  }
87 }
88 
89 
90 template<class Type>
91 void Foam::probes::sampleAndWrite
92 (
94 )
95 {
96  Field<Type> values(sample(sField));
97 
98  if (Pstream::master())
99  {
100  unsigned int w = IOstream::defaultPrecision() + 7;
101  OFstream& os = *probeFilePtrs_[sField.name()];
102 
103  os << setw(w) << sField.time().timeToUserTime(sField.time().value());
104 
105  forAll(values, probei)
106  {
107  os << ' ' << setw(w) << values[probei];
108  }
109  os << endl;
110  }
111 }
112 
113 
114 template<class Type>
115 void Foam::probes::sampleAndWrite(const fieldGroup<Type>& fields)
116 {
117  forAll(fields, fieldi)
118  {
119  if (loadFromFiles_)
120  {
121  sampleAndWrite
122  (
124  (
125  IOobject
126  (
127  fields[fieldi],
128  mesh_.time().timeName(),
129  mesh_,
132  false
133  ),
134  mesh_
135  )
136  );
137  }
138  else
139  {
140  objectRegistry::const_iterator iter = mesh_.find(fields[fieldi]);
141 
142  if
143  (
144  iter != objectRegistry::end()
145  && iter()->type()
147  )
148  {
149  sampleAndWrite
150  (
151  mesh_.lookupObject
153  (
154  fields[fieldi]
155  )
156  );
157  }
158  }
159  }
160 }
161 
162 
163 template<class Type>
164 void Foam::probes::sampleAndWriteSurfaceFields(const fieldGroup<Type>& fields)
165 {
166  forAll(fields, fieldi)
167  {
168  if (loadFromFiles_)
169  {
170  sampleAndWrite
171  (
173  (
174  IOobject
175  (
176  fields[fieldi],
177  mesh_.time().timeName(),
178  mesh_,
181  false
182  ),
183  mesh_
184  )
185  );
186  }
187  else
188  {
189  objectRegistry::const_iterator iter = mesh_.find(fields[fieldi]);
190 
191  if
192  (
193  iter != objectRegistry::end()
194  && iter()->type()
196  )
197  {
198  sampleAndWrite
199  (
200  mesh_.lookupObject
202  (
203  fields[fieldi]
204  )
205  );
206  }
207  }
208  }
209 }
210 
211 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
212 
213 template<class Type>
216 (
218 ) const
219 {
220  const Type unsetVal(-vGreat*pTraits<Type>::one);
221 
222  tmp<Field<Type>> tValues
223  (
224  new Field<Type>(this->size(), unsetVal)
225  );
226 
227  Field<Type>& values = tValues.ref();
228 
229  if (fixedLocations_)
230  {
231  autoPtr<interpolation<Type>> interpolator
232  (
233  interpolation<Type>::New(interpolationScheme_, vField)
234  );
235 
236  forAll(*this, probei)
237  {
238  if (elementList_[probei] >= 0)
239  {
240  const vector& position = operator[](probei);
241 
242  values[probei] = interpolator().interpolate
243  (
244  position,
245  elementList_[probei],
246  -1
247  );
248  }
249  }
250  }
251  else
252  {
253  forAll(*this, probei)
254  {
255  if (elementList_[probei] >= 0)
256  {
257  values[probei] = vField[elementList_[probei]];
258  }
259  }
260  }
261 
264 
265  return tValues;
266 }
267 
268 
269 template<class Type>
271 Foam::probes::sample(const word& fieldName) const
272 {
273  return sample
274  (
276  (
277  fieldName
278  )
279  );
280 }
281 
282 
283 template<class Type>
286 (
288 ) const
289 {
290  const Type unsetVal(-vGreat*pTraits<Type>::one);
291 
292  tmp<Field<Type>> tValues
293  (
294  new Field<Type>(this->size(), unsetVal)
295  );
296 
297  Field<Type>& values = tValues.ref();
298 
299  forAll(*this, probei)
300  {
301  if (faceList_[probei] >= 0)
302  {
303  values[probei] = sField[faceList_[probei]];
304  }
305  }
306 
309 
310  return tValues;
311 }
312 
313 
314 template<class Type>
317 {
318  return sample
319  (
321  (
322  fieldName
323  )
324  );
325 }
326 
327 // ************************************************************************* //
Foam::surfaceFields.
void operator()(T &x, const T &y) const
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
const word & name() const
Return name.
Definition: IOobject.H:297
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:110
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
Output to file stream.
Definition: OFstream.H:82
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Traits class for primitives.
Definition: pTraits.H:50
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:421
Generic GeometricField class.
tmp< Field< Type > > sampleSurfaceFields(const word &fieldName) const
Sample a single scalar field on all sample locations.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
tmp< Field< Type > > sample(const GeometricField< Type, fvPatchField, volMesh > &) const
Sample a volume field at all locations.
scalar y
Pre-declare SubField and related Field type.
Definition: Field.H:57
A class for handling words, derived from string.
Definition: word.H:59
friend class const_iterator
Declare friendship with the const_iterator.
Definition: HashTable.H:197
Istream and Ostream manipulators taking arguments.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
const Time & time() const
Return time.
Definition: IOobject.C:367
Abstract base class for interpolation.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.