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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "patchProbes.H"
29 #include "volFields.H"
30 #include "IOmanip.H"
31 
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 template<class Type>
37 (
39 )
40 {
41  Field<Type> values(sample(vField));
42 
43  if (Pstream::master())
44  {
45  unsigned int w = IOstream::defaultPrecision() + 7;
46  OFstream& probeStream = *probeFilePtrs_[vField.name()];
47 
48  probeStream
49  << setw(w)
50  << vField.time().timeOutputValue();
51 
52  forAll(values, probei)
53  {
54  probeStream << ' ' << setw(w) << values[probei];
55  }
56  probeStream << endl;
57  }
58 }
59 
60 
61 template<class Type>
63 (
65 )
66 {
67  Field<Type> values(sample(sField));
68 
69  if (Pstream::master())
70  {
71  unsigned int w = IOstream::defaultPrecision() + 7;
72  OFstream& probeStream = *probeFilePtrs_[sField.name()];
73 
74  probeStream
75  << setw(w)
76  << sField.time().timeOutputValue();
77 
78  forAll(values, probei)
79  {
80  probeStream << ' ' << setw(w) << values[probei];
81  }
82  probeStream << endl;
83  }
84 }
85 
86 
87 template<class Type>
89 (
91 )
92 {
93  forAll(fields, fieldi)
94  {
95  if (loadFromFiles_)
96  {
97  sampleAndWrite
98  (
100  (
101  IOobject
102  (
103  fields[fieldi],
104  mesh_.time().timeName(),
105  mesh_,
106  IOobject::MUST_READ,
107  IOobject::NO_WRITE,
108  false
109  ),
110  mesh_
111  )
112  );
113  }
114  else
115  {
116  objectRegistry::const_iterator iter = mesh_.find(fields[fieldi]);
117 
118  if
119  (
120  iter.found()
121  && iter()->type()
123  )
124  {
125  sampleAndWrite
126  (
127  mesh_.lookupObject
129  (
130  fields[fieldi]
131  )
132  );
133  }
134  }
135  }
136 }
137 
138 
139 template<class Type>
141 (
142  const fieldGroup<Type>& fields
143 )
144 {
145  forAll(fields, fieldi)
146  {
147  if (loadFromFiles_)
148  {
149  sampleAndWrite
150  (
152  (
153  IOobject
154  (
155  fields[fieldi],
156  mesh_.time().timeName(),
157  mesh_,
158  IOobject::MUST_READ,
159  IOobject::NO_WRITE,
160  false
161  ),
162  mesh_
163  )
164  );
165  }
166  else
167  {
168  objectRegistry::const_iterator iter = mesh_.find(fields[fieldi]);
169 
170  if
171  (
172  iter.found()
173  && iter()->type()
175  )
176  {
177  sampleAndWrite
178  (
179  mesh_.lookupObject
181  (
182  fields[fieldi]
183  )
184  );
185  }
186  }
187  }
188 }
189 
190 
191 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192 
193 template<class Type>
196 (
198 ) const
199 {
200  const Type unsetVal(-VGREAT*pTraits<Type>::one);
201 
202  tmp<Field<Type>> tValues
203  (
204  new Field<Type>(this->size(), unsetVal)
205  );
206 
207  Field<Type>& values = tValues.ref();
208 
209  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
210 
211  forAll(*this, probei)
212  {
213  label facei = faceList_[probei];
214 
215  if (facei >= 0)
216  {
217  label patchi = patches.whichPatch(facei);
218  label localFacei = patches[patchi].whichFace(facei);
219  values[probei] = vField.boundaryField()[patchi][localFacei];
220  }
221  }
222 
223  Pstream::listCombineGather(values, isNotEqOp<Type>());
224  Pstream::listCombineScatter(values);
225 
226  return tValues;
227 }
228 
229 
230 template<class Type>
232 Foam::patchProbes::sample(const word& fieldName) const
233 {
234  return sample
235  (
237  (
238  fieldName
239  )
240  );
241 }
242 
243 
244 template<class Type>
247 (
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  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
261 
262  forAll(*this, probei)
263  {
264  label facei = faceList_[probei];
265 
266  if (facei >= 0)
267  {
268  label patchi = patches.whichPatch(facei);
269  label localFacei = patches[patchi].whichFace(facei);
270  values[probei] = sField.boundaryField()[patchi][localFacei];
271  }
272  }
273 
276 
277  return tValues;
278 }
279 
280 
281 // ************************************************************************* //
volFields.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
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::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::probes::fieldGroup
Class used for grouping field types.
Definition: probes.H:121
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::probes::mesh_
const fvMesh & mesh_
Const reference to fvMesh.
Definition: probes.H:137
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field
Generic templated field type.
Definition: Field.H:63
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::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:812
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::patchProbes::sampleAndWriteSurfaceFields
void sampleAndWriteSurfaceFields(const fieldGroup< Type > &)
Sample and write all the surface fields of the given type.
Definition: patchProbesTemplates.C:141
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::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
patchProbes.H
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
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::patchProbes::sampleAndWrite
void sampleAndWrite(const GeometricField< Type, fvPatchField, volMesh > &)
Sample and write a particular volume field.
Definition: patchProbesTemplates.C:37
Foam::GeometricField< Type, fvPatchField, volMesh >
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
sample
Minimal example by using system/controlDict.functions: