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 Copyright (C) 2021-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "patchProbes.H"
30#include "volFields.H"
31#include "surfaceFields.H"
32#include "IOmanip.H"
33
34// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35
36template<class Type>
37void Foam::patchProbes::writeValues
38(
39 const word& fieldName,
40 const Field<Type>& values,
41 const scalar timeValue
42)
43{
44 if (Pstream::master())
45 {
46 const unsigned int w = IOstream::defaultPrecision() + 7;
47 OFstream& os = *probeFilePtrs_[fieldName];
48
49 os << setw(w) << timeValue;
50
51 for (const auto& v : values)
52 {
53 os << ' ' << setw(w) << v;
54 }
55 os << endl;
56 }
57}
58
59
60template<class GeoField>
61void Foam::patchProbes::performAction
62(
63 const fieldGroup<GeoField>& fieldNames,
64 unsigned request
65)
66{
67 for (const word& fieldName : fieldNames)
68 {
69 tmp<GeoField> tfield = getOrLoadField<GeoField>(fieldName);
70
71 if (tfield)
72 {
73 const auto& field = tfield();
74 const scalar timeValue = field.time().timeOutputValue();
75
76 Field<typename GeoField::value_type> values(sample(field));
77
78 this->storeResults(fieldName, values);
79 if (request & ACTION_WRITE)
80 {
81 this->writeValues(fieldName, values, timeValue);
82 }
83 }
84 }
85}
86
87
88// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89
90template<class Type>
93{
94 const Type unsetVal(-VGREAT*pTraits<Type>::one);
95
96 auto tvalues = tmp<Field<Type>>::New(Field<Type>(this->size(), unsetVal));
97 auto& values = tvalues.ref();
98
99 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
100 const auto& bField = vField.boundaryField();
101
102 forAll(*this, probei)
103 {
104 label facei = faceList_[probei];
105
106 if (facei >= 0)
107 {
108 label patchi = patches.whichPatch(facei);
109 label localFacei = patches[patchi].whichFace(facei);
110 values[probei] = bField[patchi][localFacei];
111 }
112 }
113
115
116 return tvalues;
117}
118
119
120template<class Type>
123{
124 const Type unsetVal(-VGREAT*pTraits<Type>::one);
125
126 auto tvalues = tmp<Field<Type>>::New(Field<Type>(this->size(), unsetVal));
127 auto& values = tvalues.ref();
128
129 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
130 const auto& bField = sField.boundaryField();
131
132 forAll(*this, probei)
133 {
134 label facei = faceList_[probei];
135
136 if (facei >= 0)
137 {
138 label patchi = patches.whichPatch(facei);
139 label localFacei = patches[patchi].whichFace(facei);
140 values[probei] = bField[patchi][localFacei];
141 }
142 }
143
145
146 return tvalues;
147}
148
149
150template<class Type>
152Foam::patchProbes::sample(const word& fieldName) const
153{
154 return sample(mesh_.lookupObject<VolumeField<Type>>(fieldName));
155}
156
157
158template<class Type>
161{
162 return sample(mesh_.lookupObject<SurfaceField<Type>>(fieldName));
163}
164
165
166// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
Generic templated field type.
Definition: Field.H:82
Minimal example by using system/controlDict.functions:
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:342
static void listCombineAllGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
tmp< Field< Type > > sampleSurfaceField(const word &fieldName) const
Sample a surface field at all locations.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
HashPtrTable< OFstream > probeFilePtrs_
Current open files (non-empty on master only)
Definition: probes.H:245
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
rDeltaTY field()
const polyBoundaryMesh & patches
OBJstream os(runTime.globalPath()/outputName)
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.