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 | Copyright (C) 2011-2015 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 //- Comparison operator for probes class
38 template<class T>
39 class isNotEqOp
40 {
41 public:
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:chould 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 }
63 
64 
65 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
66 
67 template<class Type>
69 (
71 )
72 {
73  Field<Type> values(sample(vField));
74 
75  if (Pstream::master())
76  {
77  unsigned int w = IOstream::defaultPrecision() + 7;
78  OFstream& os = *probeFilePtrs_[vField.name()];
79 
80  os << setw(w) << vField.time().timeToUserTime(vField.time().value());
81 
82  forAll(values, probeI)
83  {
84  os << ' ' << setw(w) << values[probeI];
85  }
86  os << endl;
87  }
88 }
89 
90 
91 template<class Type>
93 (
95 )
96 {
97  Field<Type> values(sample(sField));
98 
99  if (Pstream::master())
100  {
101  unsigned int w = IOstream::defaultPrecision() + 7;
102  OFstream& os = *probeFilePtrs_[sField.name()];
103 
104  os << setw(w) << sField.time().timeToUserTime(sField.time().value());
105 
106  forAll(values, probeI)
107  {
108  os << ' ' << setw(w) << values[probeI];
109  }
110  os << endl;
111  }
112 }
113 
114 
115 template<class Type>
117 {
118  forAll(fields, fieldI)
119  {
120  if (loadFromFiles_)
121  {
123  (
125  (
126  IOobject
127  (
128  fields[fieldI],
129  mesh_.time().timeName(),
130  mesh_,
133  false
134  ),
135  mesh_
136  )
137  );
138  }
139  else
140  {
142 
143  if
144  (
145  iter != objectRegistry::end()
146  && iter()->type()
148  )
149  {
151  (
154  (
155  fields[fieldI]
156  )
157  );
158  }
159  }
160  }
161 }
162 
163 
164 template<class Type>
166 {
167  forAll(fields, fieldI)
168  {
169  if (loadFromFiles_)
170  {
171  sampleAndWrite
172  (
174  (
175  IOobject
176  (
177  fields[fieldI],
178  mesh_.time().timeName(),
179  mesh_,
182  false
183  ),
184  mesh_
185  )
186  );
187  }
188  else
189  {
190  objectRegistry::const_iterator iter = mesh_.find(fields[fieldI]);
191 
192  if
193  (
194  iter != objectRegistry::end()
195  && iter()->type()
197  )
198  {
199  sampleAndWrite
200  (
201  mesh_.lookupObject
203  (
204  fields[fieldI]
205  )
206  );
207  }
208  }
209  }
210 }
211 
212 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
213 
214 template<class Type>
217 (
219 ) const
220 {
221  const Type unsetVal(-VGREAT*pTraits<Type>::one);
222 
223  tmp<Field<Type> > tValues
224  (
225  new Field<Type>(this->size(), unsetVal)
226  );
227 
228  Field<Type>& values = tValues();
229 
230  if (fixedLocations_)
231  {
232  autoPtr<interpolation<Type> > interpolator
233  (
234  interpolation<Type>::New(interpolationScheme_, vField)
235  );
236 
237  forAll(*this, probeI)
238  {
239  if (elementList_[probeI] >= 0)
240  {
241  const vector& position = operator[](probeI);
242 
243  values[probeI] = interpolator().interpolate
244  (
245  position,
246  elementList_[probeI],
247  -1
248  );
249  }
250  }
251  }
252  else
253  {
254  forAll(*this, probeI)
255  {
256  if (elementList_[probeI] >= 0)
257  {
258  values[probeI] = vField[elementList_[probeI]];
259  }
260  }
261  }
262 
265 
266  return tValues;
267 }
268 
269 
270 template<class Type>
272 Foam::probes::sample(const word& fieldName) const
273 {
274  return sample
275  (
277  (
278  fieldName
279  )
280  );
281 }
282 
283 
284 template<class Type>
287 (
289 ) const
290 {
291  const Type unsetVal(-VGREAT*pTraits<Type>::one);
292 
293  tmp<Field<Type> > tValues
294  (
295  new Field<Type>(this->size(), unsetVal)
296  );
297 
298  Field<Type>& values = tValues();
299 
300  forAll(*this, probeI)
301  {
302  if (faceList_[probeI] >= 0)
303  {
304  values[probeI] = sField[faceList_[probeI]];
305  }
306  }
307 
310 
311  return tValues;
312 }
313 
314 
315 template<class Type>
318 {
319  return sample
320  (
322  (
323  fieldName
324  )
325  );
326 }
327 
328 // ************************************************************************* //
volFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
fields
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:127
Foam::isNotEqOp
Comparison operator for probes class.
Definition: probesTemplates.C:39
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
interpolation.H
Foam::isNotEqOp::operator()
void operator()(T &x, const T &y) const
Definition: probesTemplates.C:43
probes.H
Foam::probes::fieldGroup
Class used for grouping field types.
Definition: probes.H:112
Foam::HashTable< regIOobject * >::const_iterator
friend class const_iterator
Declare friendship with the const_iterator.
Definition: HashTable.H:192
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::probes::mesh_
const fvMesh & mesh_
Const reference to fvMesh.
Definition: probes.H:132
Foam::probes::sample
tmp< Field< Type > > sample(const GeometricField< Type, fvPatchField, volMesh > &) const
Sample a volume field at all locations.
Foam::Field< Type >
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::interpolation
Abstract base class for interpolation.
Definition: mappedPatchFieldBase.H:57
Foam::probes::sampleSurfaceFields
tmp< Field< Type > > sampleSurfaceFields(const word &fieldName) const
Sample a single scalar field on all sample locations.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::probes::loadFromFiles_
const bool loadFromFiles_
Load fields from files (not from objectRegistry)
Definition: probes.H:135
Foam::HashTable::find
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Foam::OFstream
Output to file stream.
Definition: OFstream.H:81
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:282
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::probes::sampleAndWrite
void sampleAndWrite(const GeometricField< Type, fvPatchField, volMesh > &)
Sample and write a particular volume field.
Definition: probesTemplates.C:69
Foam::Vector< scalar >
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::probes::sampleAndWriteSurfaceFields
void sampleAndWriteSurfaceFields(const fieldGroup< Type > &)
Sample and write all the surface fields of the given type.
Definition: probesTemplates.C:165
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
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:417
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
y
scalar y
Definition: LISASMDCalcMethod1.H:14