fieldVisualisationBase.C 20.5 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2015-2019 OpenCFD Ltd.
6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29
     \\/     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/>.

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

// OpenFOAM includes
#include "fieldVisualisationBase.H"
#include "runTimePostProcessing.H"

30 31 32
#include "doubleVector.H"
#include "foamVtkTools.H"

33 34
// VTK includes
#include "vtkArrowSource.h"
35
#include "vtkCellDataToPointData.h"
36
#include "vtkCellData.h"
37
#include "vtkColorTransferFunction.h"
38 39 40
#include "vtkCompositeDataSet.h"
#include "vtkDataObjectTreeIterator.h"
#include "vtkFieldData.h"
41 42 43 44 45 46 47 48 49
#include "vtkGlyph3D.h"
#include "vtkLookupTable.h"
#include "vtkPointData.h"
#include "vtkPolyData.h"
#include "vtkPolyDataMapper.h"
#include "vtkRenderer.h"
#include "vtkSmartPointer.h"
#include "vtkSphereSource.h"

50
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
51

52 53
const Foam::Enum
<
54 55
    Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
    colourByType
56
>
57 58
Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
colourByTypeNames
Mark Olesen's avatar
Mark Olesen committed
59
({
60 61
    { colourByType::cbColour, "colour" },
    { colourByType::cbField, "field" },
Mark Olesen's avatar
Mark Olesen committed
62
});
63 64 65

const Foam::Enum
<
66 67
    Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
    colourMapType
68
>
69 70
Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
colourMapTypeNames
Mark Olesen's avatar
Mark Olesen committed
71
({
72 73 74
    { colourMapType::cmCoolToWarm, "coolToWarm" },
    { colourMapType::cmCoolToWarm, "blueWhiteRed" },
    { colourMapType::cmColdAndHot, "coldAndHot" },
75
    { colourMapType::cmFire, "fire" },
76
    { colourMapType::cmRainbow, "rainbow" },
77
    { colourMapType::cmGreyscale, "greyscale" },
78 79
    { colourMapType::cmGreyscale, "grayscale" },
    { colourMapType::cmXray, "xray" },
Mark Olesen's avatar
Mark Olesen committed
80
});
81 82


83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 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 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 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 208 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 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 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 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388
// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //

Foam::functionObjects::runTimePostPro::fieldVisualisationBase::fieldSummary
Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
queryFieldSummary
(
    const word& fieldName,
    vtkDataSet* dataset
)
{
    fieldSummary queried;

    if (dataset)
    {
        vtkDataArray* array;

        array = vtkDataArray::SafeDownCast
        (
            dataset->GetCellData()->GetAbstractArray(fieldName.c_str())
        );

        if (array)
        {
            queried.nComponents_ = array->GetNumberOfComponents();
            queried.association_ |= FieldAssociation::CELL_DATA;
            queried.range_ += vtk::Tools::rangeOf(array);
        }

        array = vtkDataArray::SafeDownCast
        (
            dataset->GetPointData()->GetAbstractArray(fieldName.c_str())
        );

        if (array)
        {
            queried.nComponents_ = array->GetNumberOfComponents();
            queried.association_ |= FieldAssociation::POINT_DATA;
            queried.range_ += vtk::Tools::rangeOf(array);
        }
    }

    return queried;
}


Foam::functionObjects::runTimePostPro::fieldVisualisationBase::fieldSummary
Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
queryFieldSummary
(
    const word& fieldName,
    vtkCompositeDataSet* data
)
{
    fieldSummary queried;

    auto iter = vtkSmartPointer<vtkDataObjectTreeIterator>::New();

    iter->SetDataSet(data);
    iter->VisitOnlyLeavesOn();
    iter->SkipEmptyNodesOn();

    for
    (
        iter->InitTraversal();
        !iter->IsDoneWithTraversal();
        iter->GoToNextItem()
    )
    {
        vtkDataSet* dataset = vtkDataSet::SafeDownCast
        (
            iter->GetCurrentDataObject()
        );

        if (dataset)
        {
            fieldSummary local(queryFieldSummary(fieldName, dataset));

            if (!queried.nComponents_)
            {
                queried.nComponents_ = local.nComponents_;
            }

            queried.association_ |= local.association_;
            queried.range_ += local.range_;
        }
    }

    return queried;
}


Foam::functionObjects::runTimePostPro::fieldVisualisationBase::FieldAssociation
Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
queryFieldAssociation
(
    const word& fieldName,
    vtkDataSet* dataset
)
{
    unsigned where(FieldAssociation::NO_DATA);

    if (dataset)
    {
        if (dataset->GetCellData()->HasArray(fieldName.c_str()))
        {
            where |= FieldAssociation::CELL_DATA;
        }
        if (dataset->GetPointData()->HasArray(fieldName.c_str()))
        {
            where |= FieldAssociation::POINT_DATA;
        }
    }

    return FieldAssociation(where);
}


Foam::functionObjects::runTimePostPro::fieldVisualisationBase::FieldAssociation
Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
queryFieldAssociation
(
    const word& fieldName,
    vtkCompositeDataSet* data
)
{
    unsigned where(FieldAssociation::NO_DATA);

    auto iter = vtkSmartPointer<vtkDataObjectTreeIterator>::New();

    iter->SetDataSet(data);
    iter->VisitOnlyLeavesOn();
    iter->SkipEmptyNodesOn();

    for
    (
        iter->InitTraversal();
        !iter->IsDoneWithTraversal();
        iter->GoToNextItem()
    )
    {
        vtkDataSet* dataset = vtkDataSet::SafeDownCast
        (
            iter->GetCurrentDataObject()
        );

        where |= queryFieldAssociation(fieldName, dataset);
    }

    return FieldAssociation(where);
}


void Foam::functionObjects::runTimePostPro::fieldVisualisationBase::addMagField
(
    const word& fieldName,
    vtkFieldData* fieldData
)
{
    if (!fieldData)
    {
        return;
    }

    vtkDataArray* input = vtkDataArray::SafeDownCast
    (
        fieldData->GetAbstractArray(fieldName.c_str())
    );

    if (!input)
    {
        return;
    }

    const word magFieldName = "mag(" + fieldName + ")";

    vtkDataArray* output = vtkDataArray::SafeDownCast
    (
        fieldData->GetAbstractArray(magFieldName.c_str())
    );

    if (output)
    {
        return;
    }


    // Simplfy and only handle scalar/vector input

    const int nCmpt = input->GetNumberOfComponents();
    const vtkIdType len = input->GetNumberOfTuples();

    if (nCmpt == 1)
    {
        auto data = vtkSmartPointer<vtkFloatArray>::New();

        data->SetName(magFieldName.c_str());
        data->SetNumberOfComponents(1);
        data->SetNumberOfTuples(len);

        double scratch;
        for (vtkIdType i=0; i < len; ++i)
        {
            input->GetTuple(i, &scratch);

            scratch = Foam::mag(scratch);
            data->SetTuple(i, &scratch);
        }

        fieldData->AddArray(data);
    }
    else if (nCmpt == 3)
    {
        auto data = vtkSmartPointer<vtkFloatArray>::New();

        data->SetName(magFieldName.c_str());
        data->SetNumberOfComponents(1);
        data->SetNumberOfTuples(len);

        doubleVector scratch;
        for (vtkIdType i=0; i < len; ++i)
        {
            input->GetTuple(i, scratch.v_);

            scratch.x() = Foam::mag(scratch);

            data->SetTuple(i, scratch.v_);
        }

        fieldData->AddArray(data);
    }
}


void Foam::functionObjects::runTimePostPro::fieldVisualisationBase::addMagField
(
    const word& fieldName,
    vtkDataSet* dataset
)
{
    if (dataset)
    {
        addMagField(fieldName, dataset->GetCellData());
        addMagField(fieldName, dataset->GetPointData());
    }
}


void Foam::functionObjects::runTimePostPro::fieldVisualisationBase::addMagField
(
    const word& fieldName,
    vtkCompositeDataSet* data
)
{
    auto iter = vtkSmartPointer<vtkDataObjectTreeIterator>::New();

    iter->SetDataSet(data);
    iter->VisitOnlyLeavesOn();
    iter->SkipEmptyNodesOn();

    for
    (
        iter->InitTraversal();
        !iter->IsDoneWithTraversal();
        iter->GoToNextItem()
    )
    {
        vtkDataSet* dataset = vtkDataSet::SafeDownCast
        (
            iter->GetCurrentDataObject()
        );
        addMagField(fieldName, dataset);
    }
}


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

void Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
fieldSummary::reduce()
{
    if (Pstream::parRun())
    {
        Foam::reduce(nComponents_, maxOp<int>());
        Foam::reduce(association_, bitOrOp<unsigned>());
        Foam::reduce(range_, minMaxOp<scalar>());
    }
}


Foam::Ostream& Foam::operator<<
(
    Ostream& os,
    const InfoProxy
    <
        functionObjects::runTimePostPro::fieldVisualisationBase::fieldSummary
    >& proxy
)
{
    os  << "nComponents:" << proxy.t_.nComponents_
        << " association:" << label(proxy.t_.association_)
        << " min/max:" << proxy.t_.range_;

    return os;
}


389 390
// * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //

391 392
void Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
setColourMap
393 394 395
(
    vtkLookupTable* lut
) const
396
{
397
    constexpr label nColours = 256;
398 399 400

    lut->SetNumberOfColors(nColours);

401
    auto ctf = vtkSmartPointer<vtkColorTransferFunction>::New();
402 403 404

    switch (colourMap_)
    {
405
        case cmCoolToWarm:  // ParaView: "Cool To Warm"
406 407
        {
            ctf->SetColorSpaceToDiverging();
408
            ctf->AddRGBPoint(0.0, 0.231372, 0.298039, 0.752941);
409 410
            ctf->AddRGBPoint(0.5, 0.865003, 0.865003, 0.865003);
            ctf->AddRGBPoint(1.0, 0.705882, 0.0156863, 0.14902);
411 412 413 414 415 416 417 418 419 420 421 422
            // ctf->SetNanColor(1, 1, 0);
            break;
        }

        case cmColdAndHot:  // ParaView : "Cold and Hot"
        {
            ctf->SetColorSpaceToRGB();
            ctf->AddRGBPoint(0, 0, 1, 1);
            ctf->AddRGBPoint(0.45, 0, 0, 1);
            ctf->AddRGBPoint(0.5, 0, 0, 0.5019608);
            ctf->AddRGBPoint(0.55, 1, 0, 0);
            ctf->AddRGBPoint(1, 1, 1, 0);
423 424
            break;
        }
425 426

        case cmFire:  // ParaView: Black-Body Radiation
427 428 429 430 431 432
        {
            ctf->SetColorSpaceToRGB();
            ctf->AddRGBPoint(0, 0, 0, 0);
            ctf->AddRGBPoint(0.4, 0.901961, 0, 0);
            ctf->AddRGBPoint(0.8, 0.901961, 0.901961, 0);
            ctf->AddRGBPoint(1, 1, 1, 1);
433 434 435 436 437 438 439 440 441 442 443
            // ctf->SetNanColor(0, 0.49804, 1);
            break;
        }

        case cmRainbow:
        {
            ctf->SetColorSpaceToHSV();
            ctf->AddRGBPoint(0, 0, 0, 1);
            ctf->AddRGBPoint(0.5, 0, 1, 0);
            ctf->AddRGBPoint(1, 1, 0, 0);
            // ctf->SetNanColor(0.498039, 0.498039, 0.498039);
444 445
            break;
        }
446 447

        case cmGreyscale: // ParaView: grayscale
448 449 450 451
        {
            ctf->SetColorSpaceToRGB();
            ctf->AddRGBPoint(0, 0, 0, 0);
            ctf->AddRGBPoint(1, 1, 1, 1);
452 453 454 455 456 457 458 459 460 461
            // ctf->SetNanColor(1, 0, 0);
            break;
        }

        case cmXray: // ParaView: "X ray"
        {
            ctf->SetColorSpaceToRGB();
            ctf->AddRGBPoint(0, 1, 1, 1);
            ctf->AddRGBPoint(1, 0, 0, 0);
            // ctf->SetNanColor(1, 0, 0);
462 463 464 465 466
            break;
        }
    }


467 468
    double rgba[4] = { 0, 0, 0, 1 };
    for (label i = 0; i < nColours; ++i)
469
    {
470 471
        ctf->GetColor(scalar(i)/scalar(nColours), rgba);
        lut->SetTableValue(i, rgba);
472 473 474 475
    }
}


476 477
void Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
addScalarBar
478 479 480 481 482 483
(
    const scalar position,
    vtkRenderer* renderer,
    vtkLookupTable* lut
) const
{
484 485
    // Add the scalar bar - only once!
    if (renderer && Pstream::master())
486
    {
487
        scalarBar_.add(colours_["text"]->value(position), renderer, lut);
488 489 490 491
    }
}


492 493
void Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
setField
494 495 496
(
    const scalar position,
    const word& colourFieldName,
497 498 499
    const FieldAssociation fieldAssociation,
    vtkMapper* mapper,
    vtkRenderer* renderer
500 501 502 503 504 505 506 507 508 509 510
) const
{
    mapper->InterpolateScalarsBeforeMappingOn();

    switch (colourBy_)
    {
        case cbColour:
        {
            mapper->ScalarVisibilityOff();
            break;
        }
511

512 513
        case cbField:
        {
514
            // Create look-up table for colours
515
            auto lut = vtkSmartPointer<vtkLookupTable>::New();
516 517
            setColourMap(lut);
            lut->SetVectorMode(vtkScalarsToColors::MAGNITUDE);
518
            lut->SetTableRange(range_.first(), range_.second());
519

520
            // Configure the mapper
521 522 523
            const char* fieldName = colourFieldName.c_str();
            mapper->SelectColorArray(fieldName);

524 525 526 527 528
            // Use either point or cell data
            // - if both point and cell data exists, preferentially choose
            //   point data.  This is often the case when using glyphs.

            if (fieldAssociation & FieldAssociation::POINT_DATA)
529 530 531
            {
                mapper->SetScalarModeToUsePointFieldData();
            }
532
            else if (fieldAssociation & FieldAssociation::CELL_DATA)
533 534 535 536 537 538 539 540 541 542
            {
                mapper->SetScalarModeToUseCellFieldData();
            }
            else
            {
                WarningInFunction
                    << "Unable to determine cell or point data type "
                    << "- assuming point data";
                mapper->SetScalarModeToUsePointFieldData();
            }
543 544 545 546 547
            mapper->SetScalarRange(range_.first(), range_.second());
            mapper->SetColorModeToMapScalars();
            mapper->SetLookupTable(lut);
            mapper->ScalarVisibilityOn();

548
            // Add the scalar bar
549 550 551 552 553 554 555 556 557
            addScalarBar(position, renderer, lut);
            break;
        }
    }

    mapper->Modified();
}


558 559
void Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
addGlyphs
560 561 562
(
    const scalar position,
    const word& scaleFieldName,
563
    const fieldSummary& scaleFieldInfo,
564
    const word& colourFieldName,
565
    const fieldSummary& colourFieldInfo,
566
    const scalar maxGlyphLength,
567

568 569 570 571 572
    vtkPolyData* data,
    vtkActor* actor,
    vtkRenderer* renderer
) const
{
573 574
    // Determine whether we have CellData/PointData and (scalar/vector)
    // or if we need to a cell->point data filter.
575

576
    if (!scaleFieldInfo.exists())
577 578 579 580 581 582 583
    {
        WarningInFunction
            << "Cannot add glyphs. No such cell or point field: "
            << scaleFieldName << endl;
        return;
    }

584
    if (!scaleFieldInfo.isScalar() && !scaleFieldInfo.isVector())
585 586 587 588 589 590
    {
        WarningInFunction
            << "Glyphs can only be added to scalar or vector data. "
            << "Unable to process field " << scaleFieldName << endl;
        return;
    }
591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609


    // Setup glyphs

    // The min/max data range for the input data (cell or point),
    // which will be slightly less after using a cell->point filter
    // (since it averages), but is still essentially OK.


    auto glyph = vtkSmartPointer<vtkGlyph3D>::New();
    glyph->ScalingOn();

    auto glyphMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
    glyphMapper->SetInputConnection(glyph->GetOutputPort());

    vtkSmartPointer<vtkCellDataToPointData> cellToPoint;

    // The data source is filtered or original (PointData)
    if (!scaleFieldInfo.hasPointData() || !colourFieldInfo.hasPointData())
610
    {
611 612
        // CellData - Need a cell->point filter
        cellToPoint = vtkSmartPointer<vtkCellDataToPointData>::New();
613 614
        cellToPoint->SetInputData(data);

615 616 617 618 619
        glyph->SetInputConnection(cellToPoint->GetOutputPort());
    }
    else
    {
        glyph->SetInputData(data);
620 621
    }

622

623
    if (scaleFieldInfo.nComponents_ == 1)
624
    {
625
        auto sphere = vtkSmartPointer<vtkSphereSource>::New();
626 627
        sphere->SetCenter(0, 0, 0);
        sphere->SetRadius(0.5);
628 629 630 631

        // Setting higher resolution slows the rendering significantly
        // sphere->SetPhiResolution(20);
        // sphere->SetThetaResolution(20);
632 633 634 635 636

        glyph->SetSourceConnection(sphere->GetOutputPort());

        if (maxGlyphLength > 0)
        {
637 638 639 640 641 642
            // Using range from the data:
            // glyph->SetRange
            // (
            //     scaleFieldInfo.range_.first(),
            //     scaleFieldInfo.range_.second()
            // );
643 644 645

            // Set range according to user-supplied limits
            glyph->ClampingOn();
646
            glyph->SetRange(range_.first(), range_.second());
647

648
            // If range[0] != min(value), maxGlyphLength behaviour will not
649 650 651 652 653 654 655 656 657 658 659 660
            // be correct...
            glyph->SetScaleFactor(maxGlyphLength);
        }
        else
        {
            glyph->SetScaleFactor(1);
        }
        glyph->SetScaleModeToScaleByScalar();
        glyph->OrientOff();
        glyph->SetColorModeToColorByScalar();
        glyph->SetInputArrayToProcess
        (
661 662 663
            0, // index (0) = scalars
            0, // port
            0, // connection
664
            vtkDataObject::FIELD_ASSOCIATION_POINTS,
665
            scaleFieldName.c_str()
666 667
        );
    }
668
    else if (scaleFieldInfo.nComponents_ == 3)
669
    {
670
        auto arrow = vtkSmartPointer<vtkArrowSource>::New();
671 672 673 674 675 676 677 678 679 680
        arrow->SetTipResolution(10);
        arrow->SetTipRadius(0.1);
        arrow->SetTipLength(0.35);
        arrow->SetShaftResolution(10);
        arrow->SetShaftRadius(0.03);

        glyph->SetSourceConnection(arrow->GetOutputPort());

        if (maxGlyphLength > 0)
        {
681
            // Set range according data limits
682
            glyph->ClampingOn();
683 684 685 686 687
            glyph->SetRange
            (
                scaleFieldInfo.range_.first(),
                scaleFieldInfo.range_.second()
            );
688 689 690 691 692 693 694 695 696 697 698 699
            glyph->SetScaleFactor(maxGlyphLength);
        }
        else
        {
            glyph->SetScaleFactor(1);
        }
        glyph->SetScaleModeToScaleByVector();
        glyph->OrientOn();
        glyph->SetVectorModeToUseVector();
        glyph->SetColorModeToColorByVector();
        glyph->SetInputArrayToProcess
        (
700 701 702
            1, // index (1) = vectors
            0, // port
            0, // connection
703
            vtkDataObject::FIELD_ASSOCIATION_POINTS,
704
            scaleFieldName.c_str()
705 706 707 708
        );
    }


709 710 711 712
    // Apply colouring etc.
    // We already established PointData, which as either in the original,
    // or generated with vtkCellDataToPointData filter.

713 714 715
    {
        glyph->Update();

716 717 718 719 720 721 722 723
        setField
        (
            position,
            colourFieldName,
            FieldAssociation::POINT_DATA,  // Original or after filter
            glyphMapper,
            renderer
        );
724 725 726 727 728 729 730 731 732 733 734 735

        glyphMapper->Update();

        actor->SetMapper(glyphMapper);

        renderer->AddActor(actor);
    }
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

736 737
Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
fieldVisualisationBase
738 739
(
    const dictionary& dict,
740
    const HashPtrTable<Function1<vector>>& colours
741 742 743
)
:
    colours_(colours),
744
    fieldName_(dict.get<word>("field")),
745
    smooth_(dict.getOrDefault("smooth", false)),
746 747
    colourBy_(cbColour),
    colourMap_(cmRainbow),
748 749
    range_(),
    scalarBar_()
750
{
Mark Olesen's avatar
Mark Olesen committed
751
    colourByTypeNames.readEntry("colourBy", dict, colourBy_);
752 753 754 755 756

    switch (colourBy_)
    {
        case cbColour:
        {
757
            scalarBar_.hide();
758 759
            break;
        }
760

761 762
        case cbField:
        {
763
            dict.readEntry("range", range_);
Mark Olesen's avatar
Mark Olesen committed
764
            colourMapTypeNames.readIfPresent("colourMap", dict, colourMap_);
765

766
            const dictionary* sbar = dict.findDict("scalarBar");
767

768 769 770 771 772
            if (sbar)
            {
                scalarBar_.read(*sbar);
            }
            else
773
            {
774
                scalarBar_.hide();
775
            }
776 777 778 779 780 781 782 783
            break;
        }
    }
}


// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //

784 785
Foam::functionObjects::runTimePostPro::fieldVisualisationBase::
~fieldVisualisationBase()
786 787 788 789 790
{}


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

791
const Foam::HashPtrTable<Foam::Function1<Foam::vector>, Foam::word>&
792
Foam::functionObjects::runTimePostPro::fieldVisualisationBase::colours() const
793 794 795 796 797
{
    return colours_;
}


798
const Foam::word&
799 800
Foam::functionObjects::runTimePostPro::fieldVisualisationBase::fieldName()
const
801 802 803 804 805 806
{
    return fieldName_;
}


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