patchProbesTemplates.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-2013 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 "patchProbes.H"
27 #include "volFields.H"
28 #include "IOmanip.H"
29 
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class Type>
35 (
37 )
38 {
39  Field<Type> values(sample(vField));
40 
41  if (Pstream::master())
42  {
43  unsigned int w = IOstream::defaultPrecision() + 7;
44  OFstream& probeStream = *probeFilePtrs_[vField.name()];
45 
46  probeStream
47  << setw(w)
48  << vField.time().timeToUserTime(vField.time().value());
49 
50  forAll(values, probeI)
51  {
52  probeStream << ' ' << setw(w) << values[probeI];
53  }
54  probeStream << endl;
55  }
56 }
57 
58 
59 template<class Type>
61 (
63 )
64 {
65  Field<Type> values(sample(sField));
66 
67  if (Pstream::master())
68  {
69  unsigned int w = IOstream::defaultPrecision() + 7;
70  OFstream& probeStream = *probeFilePtrs_[sField.name()];
71 
72  probeStream
73  << setw(w)
74  << sField.time().timeToUserTime(sField.time().value());
75 
76  forAll(values, probeI)
77  {
78  probeStream << ' ' << setw(w) << values[probeI];
79  }
80  probeStream << endl;
81  }
82 }
83 
84 
85 template<class Type>
87 (
89 )
90 {
91  forAll(fields, fieldI)
92  {
93  if (loadFromFiles_)
94  {
95  sampleAndWrite
96  (
98  (
99  IOobject
100  (
101  fields[fieldI],
102  mesh_.time().timeName(),
103  mesh_,
104  IOobject::MUST_READ,
105  IOobject::NO_WRITE,
106  false
107  ),
108  mesh_
109  )
110  );
111  }
112  else
113  {
114  objectRegistry::const_iterator iter = mesh_.find(fields[fieldI]);
115 
116  if
117  (
118  iter != objectRegistry::end()
119  && iter()->type()
121  )
122  {
123  sampleAndWrite
124  (
125  mesh_.lookupObject
127  (
128  fields[fieldI]
129  )
130  );
131  }
132  }
133  }
134 }
135 
136 
137 template<class Type>
139 (
140  const fieldGroup<Type>& fields
141 )
142 {
143  forAll(fields, fieldI)
144  {
145  if (loadFromFiles_)
146  {
147  sampleAndWrite
148  (
150  (
151  IOobject
152  (
153  fields[fieldI],
154  mesh_.time().timeName(),
155  mesh_,
156  IOobject::MUST_READ,
157  IOobject::NO_WRITE,
158  false
159  ),
160  mesh_
161  )
162  );
163  }
164  else
165  {
166  objectRegistry::const_iterator iter = mesh_.find(fields[fieldI]);
167 
168  if
169  (
170  iter != objectRegistry::end()
171  && iter()->type()
173  )
174  {
175  sampleAndWrite
176  (
177  mesh_.lookupObject
179  (
180  fields[fieldI]
181  )
182  );
183  }
184  }
185  }
186 }
187 
188 
189 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
190 
191 template<class Type>
194 (
196 ) const
197 {
198  const Type unsetVal(-VGREAT*pTraits<Type>::one);
199 
200  tmp<Field<Type> > tValues
201  (
202  new Field<Type>(this->size(), unsetVal)
203  );
204 
205  Field<Type>& values = tValues();
206 
207  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
208 
209  forAll(*this, probeI)
210  {
211  label faceI = faceList_[probeI];
212 
213  if (faceI >= 0)
214  {
215  label patchI = patches.whichPatch(faceI);
216  label localFaceI = patches[patchI].whichFace(faceI);
217  values[probeI] = vField.boundaryField()[patchI][localFaceI];
218  }
219  }
220 
221  Pstream::listCombineGather(values, isNotEqOp<Type>());
222  Pstream::listCombineScatter(values);
223 
224  return tValues;
225 }
226 
227 
228 template<class Type>
230 Foam::patchProbes::sample(const word& fieldName) const
231 {
232  return sample
233  (
235  (
236  fieldName
237  )
238  );
239 }
240 
241 
242 template<class Type>
245 (
247 ) const
248 {
249  const Type unsetVal(-VGREAT*pTraits<Type>::one);
250 
251  tmp<Field<Type> > tValues
252  (
253  new Field<Type>(this->size(), unsetVal)
254  );
255 
256  Field<Type>& values = tValues();
257 
258  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
259 
260  forAll(*this, probeI)
261  {
262  label faceI = faceList_[probeI];
263 
264  if (faceI >= 0)
265  {
266  label patchI = patches.whichPatch(faceI);
267  label localFaceI = patches[patchI].whichFace(faceI);
268  values[probeI] = sField.boundaryField()[patchI][localFaceI];
269  }
270  }
271 
274 
275  return tValues;
276 }
277 
278 
279 // ************************************************************************* //
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
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
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::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::probes::fieldGroup
Class used for grouping field types.
Definition: probes.H:112
Foam::probes::mesh_
const fvMesh & mesh_
Const reference to fvMesh.
Definition: probes.H:132
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field< Type >
Foam::patchProbes::sample
tmp< Field< Type > > sample(const GeometricField< Type, fvPatchField, volMesh > &) const
Sample a volume field at all locations.
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::OFstream
Output to file stream.
Definition: OFstream.H:81
Foam::patchProbes::sampleAndWriteSurfaceFields
void sampleAndWriteSurfaceFields(const fieldGroup< Type > &)
Sample and write all the surface fields of the given type.
Definition: patchProbesTemplates.C:139
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::pTraits
Traits class for primitives.
Definition: pTraits.H:50
patchProbes.H
patches
patches[0]
Definition: createSingleCellMesh.H:36
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::patchProbes::sampleAndWrite
void sampleAndWrite(const GeometricField< Type, fvPatchField, volMesh > &)
Sample and write a particular volume field.
Definition: patchProbesTemplates.C:35
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52