vtkPVFoamUpdateInfo.C 17.7 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
30
31
32
33
34
35
36

// OpenFOAM includes
#include "cellSet.H"
#include "faceSet.H"
#include "pointSet.H"
#include "IOobjectList.H"
#include "IOPtrList.H"
#include "polyBoundaryMeshEntries.H"
#include "entry.H"
#include "Cloud.H"
37
#include "vtkPVFoamReader.h"
38
39

// local headers
40
41
#include "vtkPVFoamAddToSelection.H"
#include "vtkPVFoamUpdateInfoFields.H"
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
72
73
74

// VTK includes
#include "vtkDataArraySelection.h"


// * * * * * * * * * * * * * * * Private Classes * * * * * * * * * * * * * * //

namespace Foam
{

//- A class for reading zone information without requiring a mesh
class zonesEntries
:
    public regIOobject,
    public PtrList<entry>
{

public:

    // Constructors

        explicit zonesEntries(const IOobject& io)
        :
            regIOobject(io),
            PtrList<entry>(readStream("regIOobject"))
        {
            close();
        }

   // Member functions

        bool writeData(Ostream&) const
        {
75
            NotImplemented;
76
77
78
79
80
81
82
83
84
            return false;
        }
};

}

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

template<class ZoneType>
85
Foam::wordList Foam::vtkPVFoam::getZoneNames
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
(
    const ZoneMesh<ZoneType, polyMesh>& zmesh
) const
{
    wordList names(zmesh.size());
    label nZone = 0;

    forAll(zmesh, zoneI)
    {
        if (zmesh[zoneI].size())
        {
            names[nZone++] = zmesh[zoneI].name();
        }
    }
    names.setSize(nZone);

    return names;
}


106
Foam::wordList Foam::vtkPVFoam::getZoneNames(const word& zoneType) const
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
{
    wordList names;

    // mesh not loaded - read from file
    IOobject ioObj
    (
        zoneType,
        dbPtr_().findInstance
        (
            meshDir_,
            zoneType,
            IOobject::READ_IF_PRESENT
        ),
        meshDir_,
        dbPtr_(),
        IOobject::READ_IF_PRESENT,
        IOobject::NO_WRITE,
        false
    );

127
    if (ioObj.typeHeaderOk<cellZoneMesh>(false))
128
129
130
131
132
133
134
135
136
137
138
139
140
141
    {
        zonesEntries zones(ioObj);

        names.setSize(zones.size());
        forAll(zones, zoneI)
        {
            names[zoneI] = zones[zoneI].keyword();
        }
    }

    return names;
}


142
void Foam::vtkPVFoam::updateInfoInternalMesh
143
144
145
146
147
148
(
    vtkDataArraySelection* arraySelection
)
{
    if (debug)
    {
149
        Info<< "<beg> Foam::vtkPVFoam::updateInfoInternalMesh" << endl;
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
    }

    // Determine mesh parts (internalMesh, patches...)
    //- Add internal mesh as first entry
    arrayRangeVolume_.reset(arraySelection->GetNumberOfArrays());
    arraySelection->AddArray
    (
        "internalMesh"
    );
    arrayRangeVolume_ += 1;

    if (debug)
    {
        // just for debug info
        getSelectedArrayEntries(arraySelection);

166
        Info<< "<end> Foam::vtkPVFoam::updateInfoInternalMesh" << endl;
167
168
169
170
    }
}


171
void Foam::vtkPVFoam::updateInfoLagrangian
172
173
174
175
176
177
(
    vtkDataArraySelection* arraySelection
)
{
    if (debug)
    {
178
        Info<< "<beg> Foam::vtkPVFoam::updateInfoLagrangian" << nl
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
208
209
210
211
212
213
214
215
216
            << "    " << dbPtr_->timePath()/cloud::prefix << endl;
    }


    // use the db directly since this might be called without a mesh,
    // but the region must get added back in
    fileName lagrangianPrefix(cloud::prefix);
    if (meshRegion_ != polyMesh::defaultRegion)
    {
        lagrangianPrefix = meshRegion_/cloud::prefix;
    }

    // Search for list of lagrangian objects for this time
    fileNameList cloudDirs
    (
        readDir(dbPtr_->timePath()/lagrangianPrefix, fileName::DIRECTORY)
    );

    arrayRangeLagrangian_.reset(arraySelection->GetNumberOfArrays());

    int nClouds = 0;
    forAll(cloudDirs, cloudI)
    {
        // Add cloud to GUI list
        arraySelection->AddArray
        (
            (cloudDirs[cloudI] + " - lagrangian").c_str()
        );

        ++nClouds;
    }
    arrayRangeLagrangian_ += nClouds;

    if (debug)
    {
        // just for debug info
        getSelectedArrayEntries(arraySelection);

217
        Info<< "<end> Foam::vtkPVFoam::updateInfoLagrangian" << endl;
218
219
220
221
    }
}


222
void Foam::vtkPVFoam::updateInfoPatches
223
224
225
226
227
228
229
(
    vtkDataArraySelection* arraySelection,
    stringList& enabledEntries
)
{
    if (debug)
    {
230
        Info<< "<beg> Foam::vtkPVFoam::updateInfoPatches"
231
            << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "]" << endl;
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
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
    }


    HashSet<string> enabledEntriesSet(enabledEntries);

    arrayRangePatches_.reset(arraySelection->GetNumberOfArrays());

    int nPatches = 0;
    if (meshPtr_)
    {
        const polyBoundaryMesh& patches = meshPtr_->boundaryMesh();
        const HashTable<labelList, word>& groups = patches.groupPatchIDs();
        const wordList allPatchNames = patches.names();

        // Add patch groups
        // ~~~~~~~~~~~~~~~~

        for
        (
            HashTable<labelList, word>::const_iterator iter = groups.begin();
            iter != groups.end();
            ++iter
        )
        {
            const word& groupName = iter.key();
            const labelList& patchIDs = iter();

            label nFaces = 0;
            forAll(patchIDs, i)
            {
                nFaces += patches[patchIDs[i]].size();
            }

            // Valid patch if nFace > 0 - add patch to GUI list
            if (nFaces)
            {
                string vtkGrpName = groupName + " - group";
                arraySelection->AddArray(vtkGrpName.c_str());

                ++nPatches;

                if (enabledEntriesSet.found(vtkGrpName))
                {
                    if (!reader_->GetShowGroupsOnly())
                    {
                        enabledEntriesSet.erase(vtkGrpName);
                        forAll(patchIDs, i)
                        {
                            const polyPatch& pp = patches[patchIDs[i]];
                            if (pp.size())
                            {
                                string vtkPatchName = pp.name() + " - patch";
                                enabledEntriesSet.insert(vtkPatchName);
                            }
                        }
                    }
                }
            }
        }


        // Add patches
        // ~~~~~~~~~~~

        if (!reader_->GetShowGroupsOnly())
        {
298
            forAll(patches, patchi)
299
            {
300
                const polyPatch& pp = patches[patchi];
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

                if (pp.size())
                {
                    // Add patch to GUI list
                    arraySelection->AddArray
                    (
                        (pp.name() + " - patch").c_str()
                    );

                    ++nPatches;
                }
            }
        }
    }
    else
    {
        // mesh not loaded - read from file
        // but this could fail if we've supplied a bad region name
        IOobject ioObj
        (
            "boundary",
            dbPtr_().findInstance
            (
                meshDir_,
                "boundary",
                IOobject::READ_IF_PRESENT
            ),
            meshDir_,
            dbPtr_(),
            IOobject::READ_IF_PRESENT,
            IOobject::NO_WRITE,
            false
        );

        // this should only ever fail if the mesh region doesn't exist
336
        if (ioObj.typeHeaderOk<polyBoundaryMesh>(true))
337
338
339
340
341
342
343
344
345
346
        {
            polyBoundaryMeshEntries patchEntries(ioObj);


            // Read patches and determine sizes
            // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

            wordList names(patchEntries.size());
            labelList sizes(patchEntries.size());

347
            forAll(patchEntries, patchi)
348
            {
349
                const dictionary& patchDict = patchEntries[patchi].dict();
350

351
352
                sizes[patchi] = readLabel(patchDict.lookup("nFaces"));
                names[patchi] = patchEntries[patchi].keyword();
353
354
355
356
357
358
359
360
            }


            // Add (non-zero) patch groups to the list of mesh parts
            // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

            HashTable<labelList, word> groups(patchEntries.size());

361
            forAll(patchEntries, patchi)
362
            {
363
                const dictionary& patchDict = patchEntries[patchi].dict();
364
365
366
367
368
369
370
371
372
373
374
375

                wordList groupNames;
                patchDict.readIfPresent("inGroups", groupNames);

                forAll(groupNames, groupI)
                {
                    HashTable<labelList, word>::iterator iter = groups.find
                    (
                        groupNames[groupI]
                    );
                    if (iter != groups.end())
                    {
376
                        iter().append(patchi);
377
378
379
                    }
                    else
                    {
380
                        groups.insert(groupNames[groupI], labelList(1, patchi));
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
                    }
                }
            }

            for
            (
                HashTable<labelList, word>::const_iterator iter =
                    groups.begin();
                iter != groups.end();
                ++iter
            )
            {
                const word& groupName = iter.key();
                const labelList& patchIDs = iter();

                label nFaces = 0;
                forAll(patchIDs, i)
                {
                    nFaces += sizes[patchIDs[i]];
                }

                // Valid patch if nFace > 0 - add patch to GUI list
                if (nFaces)
                {
                    string vtkGrpName = groupName + " - group";
                    arraySelection->AddArray(vtkGrpName.c_str());

                    ++nPatches;

                    if (enabledEntriesSet.found(vtkGrpName))
                    {
                        if (!reader_->GetShowGroupsOnly())
                        {
                            enabledEntriesSet.erase(vtkGrpName);
                            forAll(patchIDs, i)
                            {
                                if (sizes[patchIDs[i]])
                                {
                                    string vtkPatchName =
                                        names[patchIDs[i]] + " - patch";
                                    enabledEntriesSet.insert(vtkPatchName);
                                }
                            }
                        }
                    }
                }
            }


            // Add (non-zero) patches to the list of mesh parts
            // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

            if (!reader_->GetShowGroupsOnly())
            {
435
                forAll(names, patchi)
436
437
                {
                    // Valid patch if nFace > 0 - add patch to GUI list
438
                    if (sizes[patchi])
439
440
441
                    {
                        arraySelection->AddArray
                        (
442
                            (names[patchi] + " - patch").c_str()
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
                        );

                        ++nPatches;
                    }
                }
            }
        }
    }
    arrayRangePatches_ += nPatches;

    // Update enabled entries in case of group selection
    enabledEntries = enabledEntriesSet.toc();

    if (debug)
    {
        // just for debug info
        getSelectedArrayEntries(arraySelection);

461
        Info<< "<end> Foam::vtkPVFoam::updateInfoPatches" << endl;
462
463
464
465
    }
}


466
void Foam::vtkPVFoam::updateInfoZones
467
468
469
470
471
472
473
474
475
476
477
(
    vtkDataArraySelection* arraySelection
)
{
    if (!reader_->GetIncludeZones())
    {
        return;
    }

    if (debug)
    {
478
        Info<< "<beg> Foam::vtkPVFoam::updateInfoZones"
479
            << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "]" << endl;
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
    }

    wordList namesLst;

    //
    // cellZones information
    // ~~~~~~~~~~~~~~~~~~~~~
    if (meshPtr_)
    {
        namesLst = getZoneNames(meshPtr_->cellZones());
    }
    else
    {
        namesLst = getZoneNames("cellZones");
    }

    arrayRangeCellZones_.reset(arraySelection->GetNumberOfArrays());
    forAll(namesLst, elemI)
    {
        arraySelection->AddArray
        (
            (namesLst[elemI] + " - cellZone").c_str()
        );
    }
    arrayRangeCellZones_ += namesLst.size();


    //
    // faceZones information
    // ~~~~~~~~~~~~~~~~~~~~~
    if (meshPtr_)
    {
        namesLst = getZoneNames(meshPtr_->faceZones());
    }
    else
    {
        namesLst = getZoneNames("faceZones");
    }

    arrayRangeFaceZones_.reset(arraySelection->GetNumberOfArrays());
    forAll(namesLst, elemI)
    {
        arraySelection->AddArray
        (
            (namesLst[elemI] + " - faceZone").c_str()
        );
    }
    arrayRangeFaceZones_ += namesLst.size();


    //
    // pointZones information
    // ~~~~~~~~~~~~~~~~~~~~~~
    if (meshPtr_)
    {
        namesLst = getZoneNames(meshPtr_->pointZones());
    }
    else
    {
        namesLst = getZoneNames("pointZones");
    }

    arrayRangePointZones_.reset(arraySelection->GetNumberOfArrays());
    forAll(namesLst, elemI)
    {
        arraySelection->AddArray
        (
            (namesLst[elemI] + " - pointZone").c_str()
        );
    }
    arrayRangePointZones_ += namesLst.size();

    if (debug)
    {
        // just for debug info
        getSelectedArrayEntries(arraySelection);

557
        Info<< "<end> Foam::vtkPVFoam::updateInfoZones" << endl;
558
559
560
561
    }
}


562
void Foam::vtkPVFoam::updateInfoSets
563
564
565
566
567
568
569
570
571
572
573
(
    vtkDataArraySelection* arraySelection
)
{
    if (!reader_->GetIncludeSets())
    {
        return;
    }

    if (debug)
    {
574
        Info<< "<beg> Foam::vtkPVFoam::updateInfoSets" << endl;
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
    }

    // Add names of sets. Search for last time directory with a sets
    // subdirectory. Take care not to search beyond the last mesh.

    word facesInstance = dbPtr_().findInstance
    (
        meshDir_,
        "faces",
        IOobject::READ_IF_PRESENT
    );

    word setsInstance = dbPtr_().findInstance
    (
        meshDir_/"sets",
        word::null,
        IOobject::READ_IF_PRESENT,
        facesInstance
    );

    IOobjectList objects(dbPtr_(), setsInstance, meshDir_/"sets");

    if (debug)
    {
599
        Info<< "     Foam::vtkPVFoam::updateInfoSets read "
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
            << objects.names() << " from " << setsInstance << endl;
    }


    arrayRangeCellSets_.reset(arraySelection->GetNumberOfArrays());
    arrayRangeCellSets_ += addToSelection<cellSet>
    (
        arraySelection,
        objects,
        " - cellSet"
    );

    arrayRangeFaceSets_.reset(arraySelection->GetNumberOfArrays());
    arrayRangeFaceSets_ += addToSelection<faceSet>
    (
        arraySelection,
        objects,
        " - faceSet"
    );

    arrayRangePointSets_.reset(arraySelection->GetNumberOfArrays());
    arrayRangePointSets_ += addToSelection<pointSet>
    (
        arraySelection,
        objects,
        " - pointSet"
    );

    if (debug)
    {
        // just for debug info
        getSelectedArrayEntries(arraySelection);

633
        Info<< "<end> Foam::vtkPVFoam::updateInfoSets" << endl;
634
635
636
637
    }
}


638
void Foam::vtkPVFoam::updateInfoLagrangianFields()
639
640
641
{
    if (debug)
    {
642
        Info<< "<beg> Foam::vtkPVFoam::updateInfoLagrangianFields"
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
            << endl;
    }

    vtkDataArraySelection* fieldSelection =
        reader_->GetLagrangianFieldSelection();

    // preserve the enabled selections
    stringList enabledEntries = getSelectedArrayEntries(fieldSelection);
    fieldSelection->RemoveAllArrays();

    // TODO - currently only get fields from ONE cloud
    // have to decide if the second set of fields get mixed in
    // or dealt with separately

    const arrayRange& range = arrayRangeLagrangian_;
    if (range.empty())
    {
        return;
    }

    int partId = range.start();
    word cloudName = getPartName(partId);

    // use the db directly since this might be called without a mesh,
    // but the region must get added back in
    fileName lagrangianPrefix(cloud::prefix);
    if (meshRegion_ != polyMesh::defaultRegion)
    {
        lagrangianPrefix = meshRegion_/cloud::prefix;
    }

    IOobjectList objects
    (
        dbPtr_(),
        dbPtr_().timeName(),
        lagrangianPrefix/cloudName
    );

681
    addToSelection<IOField<label>>
682
683
684
685
    (
        fieldSelection,
        objects
    );
686
    addToSelection<IOField<scalar>>
687
688
689
690
    (
        fieldSelection,
        objects
    );
691
    addToSelection<IOField<vector>>
692
693
694
695
    (
        fieldSelection,
        objects
    );
696
    addToSelection<IOField<sphericalTensor>>
697
698
699
700
701
    (
        fieldSelection,

        objects
    );
702
    addToSelection<IOField<symmTensor>>
703
704
705
706
    (
        fieldSelection,
        objects
    );
707
    addToSelection<IOField<tensor>>
708
709
710
711
712
713
714
715
716
717
    (
        fieldSelection,
        objects
    );

    // restore the enabled selections
    setSelectedArrayEntries(fieldSelection, enabledEntries);

    if (debug)
    {
718
        Info<< "<end> Foam::vtkPVFoam::updateInfoLagrangianFields - "
719
720
721
722
723
724
            << "lagrangian objects.size() = " << objects.size() << endl;
    }
}


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