vtkPVFoamFields.C 8.14 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
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    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.

    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
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

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

26
#include "vtkPVFoam.H"
27
28
29

// OpenFOAM includes
#include "IOobjectList.H"
30
#include "vtkPVFoamReader.h"
31
32
33
34
35
36
37
38

// VTK includes
#include "vtkDataArraySelection.h"
#include "vtkPolyData.h"
#include "vtkUnstructuredGrid.h"

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

39
40
41
#include "vtkPVFoamVolFields.H"
#include "vtkPVFoamPointFields.H"
#include "vtkPVFoamLagrangianFields.H"
42
43


44
void Foam::vtkPVFoam::pruneObjectList
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
(
    IOobjectList& objects,
    const wordHashSet& selected
)
{
    // hash all the selected field names
    if (selected.empty())
    {
        objects.clear();
    }

    // only keep selected fields
    forAllIter(IOobjectList, objects, iter)
    {
        if (!selected.found(iter()->name()))
        {
            objects.erase(iter);
        }
    }
}


67
void Foam::vtkPVFoam::convertVolFields
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
(
    vtkMultiBlockDataSet* output
)
{
    const fvMesh& mesh = *meshPtr_;

    wordHashSet selectedFields = getSelected
    (
        reader_->GetVolFieldSelection()
    );

    if (selectedFields.empty())
    {
        return;
    }

    // Get objects (fields) for this time - only keep selected fields
    // the region name is already in the mesh db
    IOobjectList objects(mesh, dbPtr_().timeName());
    pruneObjectList(objects, selectedFields);

    if (objects.empty())
    {
        return;
    }

    if (debug)
    {
96
        Info<< "<beg> Foam::vtkPVFoam::convertVolFields" << nl
97
98
99
100
101
102
103
104
105
106
            << "converting OpenFOAM volume fields" << endl;
        forAllConstIter(IOobjectList, objects, iter)
        {
            Info<< "  " << iter()->name()
                << " == " << iter()->objectPath() << nl;
        }
        printMemory();
    }


107
    PtrList<PrimitivePatchInterpolation<primitivePatch>>
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
        ppInterpList(mesh.boundaryMesh().size());

    forAll(ppInterpList, i)
    {
        ppInterpList.set
        (
            i,
            new PrimitivePatchInterpolation<primitivePatch>
            (
                mesh.boundaryMesh()[i]
            )
        );
    }


    bool interpFields = reader_->GetInterpolateVolFields();

    convertVolFields<scalar>
    (
        mesh, ppInterpList, objects, interpFields, output
    );
    convertVolFields<vector>
    (
        mesh, ppInterpList, objects, interpFields, output
    );
    convertVolFields<sphericalTensor>
    (
        mesh, ppInterpList, objects, interpFields, output
    );
    convertVolFields<symmTensor>
    (
        mesh, ppInterpList, objects, interpFields, output
    );
    convertVolFields<tensor>
    (
        mesh, ppInterpList, objects, interpFields, output
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
    );

    convertDimFields<scalar>
    (
        mesh, ppInterpList, objects, interpFields, output
    );
    convertDimFields<vector>
    (
        mesh, ppInterpList, objects, interpFields, output
    );
    convertDimFields<sphericalTensor>
    (
        mesh, ppInterpList, objects, interpFields, output
    );
    convertDimFields<symmTensor>
    (
        mesh, ppInterpList, objects, interpFields, output
    );
    convertDimFields<tensor>
    (
        mesh, ppInterpList, objects, interpFields, output
165
166
167
168
    );

    if (debug)
    {
169
        Info<< "<end> Foam::vtkPVFoam::convertVolFields" << endl;
170
171
172
173
174
        printMemory();
    }
}


175
void Foam::vtkPVFoam::convertPointFields
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
(
    vtkMultiBlockDataSet* output
)
{
    const fvMesh& mesh = *meshPtr_;

    wordHashSet selectedFields = getSelected
    (
        reader_->GetPointFieldSelection()
    );

    if (selectedFields.empty())
    {
        if (debug)
        {
            Info<< "no point fields selected" << endl;
        }
        return;
    }

    // Get objects (fields) for this time - only keep selected fields
    // the region name is already in the mesh db
    IOobjectList objects(mesh, dbPtr_().timeName());
    pruneObjectList(objects, selectedFields);

    if (objects.empty())
    {
        return;
    }

    if (debug)
    {
208
        Info<< "<beg> Foam::vtkPVFoam::convertPointFields" << nl
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
            << "converting OpenFOAM volume fields -> point fields" << endl;
        forAllConstIter(IOobjectList, objects, iter)
        {
            Info<< "  " << iter()->name()
                << " == " << iter()->objectPath() << nl;
        }
        printMemory();
    }

    // Construct interpolation on the raw mesh
    const pointMesh& pMesh = pointMesh::New(mesh);


    convertPointFields<scalar>
    (
        mesh, pMesh, objects, output
    );
    convertPointFields<vector>
    (
        mesh, pMesh, objects, output
    );
    convertPointFields<sphericalTensor>
    (
        mesh, pMesh, objects, output
    );
    convertPointFields<symmTensor>
    (
        mesh, pMesh, objects, output
    );
    convertPointFields<tensor>
    (
        mesh, pMesh, objects, output
    );

    if (debug)
    {
245
        Info<< "<end> Foam::vtkPVFoam::convertPointFields" << endl;
246
247
248
249
250
        printMemory();
    }
}


251
void Foam::vtkPVFoam::convertLagrangianFields
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
(
    vtkMultiBlockDataSet* output
)
{
    arrayRange& range = arrayRangeLagrangian_;
    const fvMesh& mesh = *meshPtr_;

    wordHashSet selectedFields = getSelected
    (
        reader_->GetLagrangianFieldSelection()
    );

    if (selectedFields.empty())
    {
        return;
    }

    if (debug)
    {
271
        Info<< "<beg> Foam::vtkPVFoam::convertLagrangianFields" << endl;
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
        printMemory();
    }

    for (int partId = range.start(); partId < range.end(); ++partId)
    {
        const word  cloudName = getPartName(partId);
        const label datasetNo = partDataset_[partId];

        if (!partStatus_[partId] || datasetNo < 0)
        {
            continue;
        }


        // Get the Lagrangian fields for this time and this cloud
        // but only keep selected fields
        // the region name is already in the mesh db
        IOobjectList objects
        (
            mesh,
            dbPtr_().timeName(),
            cloud::prefix/cloudName
        );
        pruneObjectList(objects, selectedFields);

        if (objects.empty())
        {
            continue;
        }

        if (debug)
        {
            Info<< "converting OpenFOAM lagrangian fields" << nl;
            forAllConstIter(IOobjectList, objects, iter)
            {
                Info<< "  " << iter()->name()
                    << " == " << iter()->objectPath() << nl;
            }
        }

        convertLagrangianFields<label>
        (
            objects, output, datasetNo
        );
        convertLagrangianFields<scalar>
        (
            objects, output, datasetNo
        );
        convertLagrangianFields<vector>
        (
            objects, output, datasetNo
        );
        convertLagrangianFields<sphericalTensor>
        (
            objects, output, datasetNo
        );
        convertLagrangianFields<symmTensor>
        (
            objects, output, datasetNo
        );
        convertLagrangianFields<tensor>
        (
            objects, output, datasetNo
        );
    }

    if (debug)
    {
340
        Info<< "<end> Foam::vtkPVFoam::convertLagrangianFields" << endl;
341
342
343
344
345
346
        printMemory();
    }
}


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