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-2022 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 (
69  const VolField<Type>& vField
70 )
71 {
72  Field<Type> values(sample(vField));
73 
74  if (Pstream::master())
75  {
76  const unsigned int w = IOstream::defaultPrecision() + 7;
77  OFstream& os = *probeFilePtrs_[vField.name()];
78 
79  os << setw(w) << vField.time().userTimeValue();
80 
81  forAll(values, probei)
82  {
83  OStringStream buf;
84  buf << values[probei];
85  os << ' ' << setw(w) << buf.str().c_str();
86  }
87  os << endl;
88  }
89 }
90 
91 
92 template<class Type>
93 void Foam::probes::sampleAndWrite
94 (
95  const SurfaceField<Type>& sField
96 )
97 {
98  Field<Type> values(sample(sField));
99 
100  if (Pstream::master())
101  {
102  const unsigned int w = IOstream::defaultPrecision() + 7;
103  OFstream& os = *probeFilePtrs_[sField.name()];
104 
105  os << sField.time().userTimeValue();
106 
107  forAll(values, probei)
108  {
109  OStringStream buf;
110  buf << values[probei];
111  os << ' ' << setw(w) << buf.str().c_str();
112  }
113  os << endl;
114  }
115 }
116 
117 
118 template<class Type>
119 void Foam::probes::sampleAndWrite(const fieldGroup<Type>& fields)
120 {
121  forAll(fields, fieldi)
122  {
123  objectRegistry::const_iterator iter = mesh_.find(fields[fieldi]);
124 
125  if
126  (
127  iter != objectRegistry::end()
128  && iter()->type()
130  )
131  {
132  sampleAndWrite
133  (
134  mesh_.lookupObject
135  <VolField<Type>>
136  (
137  fields[fieldi]
138  )
139  );
140  }
141  }
142 }
143 
144 
145 template<class Type>
146 void Foam::probes::sampleAndWriteSurfaceFields(const fieldGroup<Type>& fields)
147 {
148  forAll(fields, fieldi)
149  {
150  objectRegistry::const_iterator iter = mesh_.find(fields[fieldi]);
151 
152  if
153  (
154  iter != objectRegistry::end()
155  && iter()->type()
157  )
158  {
159  sampleAndWrite
160  (
161  mesh_.lookupObject
162  <SurfaceField<Type>>
163  (
164  fields[fieldi]
165  )
166  );
167  }
168  }
169 }
170 
171 
172 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
173 
174 template<class Type>
177 (
178  const VolField<Type>& vField
179 ) const
180 {
181  const Type unsetVal(-vGreat*pTraits<Type>::one);
182 
183  tmp<Field<Type>> tValues
184  (
185  new Field<Type>(this->size(), unsetVal)
186  );
187 
188  Field<Type>& values = tValues.ref();
189 
190  if (fixedLocations_)
191  {
192  autoPtr<interpolation<Type>> interpolator
193  (
194  interpolation<Type>::New(interpolationScheme_, vField)
195  );
196 
197  forAll(*this, probei)
198  {
199  if (elementList_[probei] >= 0)
200  {
201  const vector& position = operator[](probei);
202 
203  values[probei] = interpolator().interpolate
204  (
205  position,
206  elementList_[probei],
207  -1
208  );
209  }
210  }
211  }
212  else
213  {
214  forAll(*this, probei)
215  {
216  if (elementList_[probei] >= 0)
217  {
218  values[probei] = vField[elementList_[probei]];
219  }
220  }
221  }
222 
225 
226  return tValues;
227 }
228 
229 
230 template<class Type>
232 Foam::probes::sample(const word& fieldName) const
233 {
234  return sample
235  (
236  mesh_.lookupObject<VolField<Type>>
237  (
238  fieldName
239  )
240  );
241 }
242 
243 
244 template<class Type>
247 (
248  const SurfaceField<Type>& sField
249 ) const
250 {
251  const Type unsetVal(-vGreat*pTraits<Type>::one);
252 
253  tmp<Field<Type>> tValues
254  (
255  new Field<Type>(this->size(), unsetVal)
256  );
257 
258  Field<Type>& values = tValues.ref();
259 
260  forAll(*this, probei)
261  {
262  if (faceList_[probei] >= 0)
263  {
264  values[probei] = sField[faceList_[probei]];
265  }
266  }
267 
270 
271  return tValues;
272 }
273 
274 
275 template<class Type>
278 {
279  return sample
280  (
281  mesh_.lookupObject<SurfaceField<Type>>
282  (
283  fieldName
284  )
285  );
286 }
287 
288 // ************************************************************************* //
Istream and Ostream manipulators taking arguments.
scalar y
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Pre-declare SubField and related Field type.
Definition: Field.H:82
static const char *const typeName
Definition: Field.H:105
Generic GeometricField class.
friend class const_iterator
Declare friendship with the const_iterator.
Definition: HashTable.H:197
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:458
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Abstract base class for interpolation.
Definition: interpolation.H:55
void operator()(T &x, const T &y) const
Traits class for primitives.
Definition: pTraits.H:53
HashPtrTable< OFstream > probeFilePtrs_
Current open files.
Definition: probes.H:132
tmp< Field< Type > > sample(const VolField< Type > &) const
Sample a volume field at all locations.
tmp< Field< Type > > sampleSurfaceFields(const word &fieldName) const
Sample a single scalar field on all sample locations.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
Foam::surfaceFields.