probesTemplates.C 7.75 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
6
7
8
9
10
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

11
12
13
14
    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.
15
16
17
18
19
20
21

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    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
27

\*---------------------------------------------------------------------------*/

#include "probes.H"
#include "volFields.H"
28
#include "surfaceFields.H"
29
#include "IOmanip.H"
30
#include "interpolation.H"
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

namespace Foam
{

template<class T>
class isNotEqOp
{
public:

    void operator()(T& x, const T& y) const
    {
        const T unsetVal(-VGREAT*pTraits<T>::one);

        if (x != unsetVal)
        {
            // Keep x.

            // Note:chould check for y != unsetVal but multiple sample cells
            // already handled in read().
        }
        else
        {
            // x is not set. y might be.
            x = y;
        }
    }
};

}


// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //

template<class Type>
void Foam::probes::sampleAndWrite
(
    const GeometricField<Type, fvPatchField, volMesh>& vField
)
{
72
    Field<Type> values(sample(vField));
73
74
75
76

    if (Pstream::master())
    {
        unsigned int w = IOstream::defaultPrecision() + 7;
77
        OFstream& os = *probeFilePtrs_[vField.name()];
78

79
        os  << setw(w) << vField.time().timeOutputValue();
80

81
        forAll(values, probei)
82
        {
83
            os  << ' ' << setw(w) << values[probei];
84
        }
85
        os  << endl;
86
87
88
89
    }
}


90
template<class Type>
91
92
void Foam::probes::sampleAndWrite
(
93
    const GeometricField<Type, fvsPatchField, surfaceMesh>& sField
94
)
95
{
96
    Field<Type> values(sample(sField));
97
98
99
100

    if (Pstream::master())
    {
        unsigned int w = IOstream::defaultPrecision() + 7;
101
        OFstream& os = *probeFilePtrs_[sField.name()];
102

103
        os  << setw(w) << sField.time().timeOutputValue();
104

105
        forAll(values, probei)
106
        {
107
            os  << ' ' << setw(w) << values[probei];
108
109
110
111
112
113
        }
        os  << endl;
    }
}


Henry's avatar
Henry committed
114
template<class Type>
115
void Foam::probes::sampleAndWrite(const fieldGroup<Type>& fields)
116
{
117
    forAll(fields, fieldi)
118
119
120
121
122
123
124
125
126
    {
        if (loadFromFiles_)
        {
            sampleAndWrite
            (
                GeometricField<Type, fvPatchField, volMesh>
                (
                    IOobject
                    (
127
                        fields[fieldi],
128
129
                        mesh_.time().timeName(),
                        mesh_,
130
131
132
133
                        IOobject::MUST_READ,
                        IOobject::NO_WRITE,
                        false
                    ),
134
                    mesh_
135
136
137
138
139
                )
            );
        }
        else
        {
140
            objectRegistry::const_iterator iter = mesh_.find(fields[fieldi]);
141
142
143

            if
            (
144
                iter.found()
145
146
147
148
149
150
             && iter()->type()
             == GeometricField<Type, fvPatchField, volMesh>::typeName
            )
            {
                sampleAndWrite
                (
151
                    mesh_.lookupObject
152
                    <GeometricField<Type, fvPatchField, volMesh>>
153
                    (
154
                        fields[fieldi]
155
156
157
158
159
160
161
162
                    )
                );
            }
        }
    }
}


163
164
165
template<class Type>
void Foam::probes::sampleAndWriteSurfaceFields(const fieldGroup<Type>& fields)
{
166
    forAll(fields, fieldi)
167
168
169
170
171
172
173
174
175
    {
        if (loadFromFiles_)
        {
            sampleAndWrite
            (
                GeometricField<Type, fvsPatchField, surfaceMesh>
                (
                    IOobject
                    (
176
                        fields[fieldi],
177
178
179
180
181
182
183
184
185
186
187
188
                        mesh_.time().timeName(),
                        mesh_,
                        IOobject::MUST_READ,
                        IOobject::NO_WRITE,
                        false
                    ),
                    mesh_
                )
            );
        }
        else
        {
189
            objectRegistry::const_iterator iter = mesh_.find(fields[fieldi]);
190
191
192

            if
            (
193
                iter.found()
194
195
196
197
198
199
200
             && iter()->type()
             == GeometricField<Type, fvsPatchField, surfaceMesh>::typeName
            )
            {
                sampleAndWrite
                (
                    mesh_.lookupObject
201
                    <GeometricField<Type, fvsPatchField, surfaceMesh>>
202
                    (
203
                        fields[fieldi]
204
205
206
207
208
209
210
                    )
                );
            }
        }
    }
}

211
212
213
// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

template<class Type>
214
Foam::tmp<Foam::Field<Type>>
215
216
217
218
219
220
221
Foam::probes::sample
(
    const GeometricField<Type, fvPatchField, volMesh>& vField
) const
{
    const Type unsetVal(-VGREAT*pTraits<Type>::one);

222
    tmp<Field<Type>> tValues
223
    (
224
        new Field<Type>(this->size(), unsetVal)
225
226
    );

227
    Field<Type>& values = tValues.ref();
228

229
    if (fixedLocations_)
230
    {
231
        autoPtr<interpolation<Type>> interpolator
232
233
234
235
        (
            interpolation<Type>::New(interpolationScheme_, vField)
        );

236
        forAll(*this, probei)
237
        {
238
            if (elementList_[probei] >= 0)
239
            {
240
                const vector& position = operator[](probei);
241

242
                values[probei] = interpolator().interpolate
243
244
                (
                    position,
245
                    elementList_[probei],
246
247
248
249
250
251
252
                    -1
                );
            }
        }
    }
    else
    {
253
        forAll(*this, probei)
254
        {
255
            if (elementList_[probei] >= 0)
256
            {
257
                values[probei] = vField[elementList_[probei]];
258
            }
259
260
261
262
263
264
265
266
267
268
269
        }
    }

    Pstream::listCombineGather(values, isNotEqOp<Type>());
    Pstream::listCombineScatter(values);

    return tValues;
}


template<class Type>
270
Foam::tmp<Foam::Field<Type>>
271
272
273
274
Foam::probes::sample(const word& fieldName) const
{
    return sample
    (
275
        mesh_.lookupObject<GeometricField<Type, fvPatchField, volMesh>>
276
277
278
279
280
281
282
        (
            fieldName
        )
    );
}


283
template<class Type>
284
Foam::tmp<Foam::Field<Type>>
285
286
Foam::probes::sample
(
287
    const GeometricField<Type, fvsPatchField, surfaceMesh>& sField
288
289
290
291
) const
{
    const Type unsetVal(-VGREAT*pTraits<Type>::one);

292
    tmp<Field<Type>> tValues
293
294
295
296
    (
        new Field<Type>(this->size(), unsetVal)
    );

297
    Field<Type>& values = tValues.ref();
298

299
    forAll(*this, probei)
300
    {
301
        if (faceList_[probei] >= 0)
302
        {
303
            values[probei] = sField[faceList_[probei]];
304
305
306
307
308
309
310
311
312
313
314
        }
    }

    Pstream::listCombineGather(values, isNotEqOp<Type>());
    Pstream::listCombineScatter(values);

    return tValues;
}


template<class Type>
315
Foam::tmp<Foam::Field<Type>>
316
317
318
319
Foam::probes::sampleSurfaceFields(const word& fieldName) const
{
    return sample
    (
320
        mesh_.lookupObject<GeometricField<Type, fvsPatchField, surfaceMesh>>
321
322
323
324
325
326
        (
            fieldName
        )
    );
}

327
// ************************************************************************* //