Skip to content
Snippets Groups Projects
runTimePostProcessing.C 6.21 KiB
Newer Older
  • Learn to ignore specific revisions
  • /*---------------------------------------------------------------------------*\
      =========                 |
      \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
       \\    /   O peration     |
    
        \\  /    A nd           | Copyright (C) 2015-2016 OpenCFD Ltd.
    
         \\/     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 "runTimePostProcessing.H"
    #include "dictionary.H"
    #include "pointData.H"
    #include "pathline.H"
    #include "surface.H"
    #include "text.H"
    #include "Time.H"
    
    #include "sigFpe.H"
    
    #include "addToRunTimeSelectionTable.H"
    
    
    // VTK includes
    #include "vtkPolyDataMapper.h"
    #include "vtkRenderer.h"
    #include "vtkRenderWindow.h"
    #include "vtkSmartPointer.h"
    
    #include "vtkLight.h"
    
    // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
    
    namespace Foam
    
    {
        defineTypeNameAndDebug(runTimePostProcessing, 0);
    
    
        addToRunTimeSelectionTable
        (
            functionObject,
            runTimePostProcessing,
            dictionary
        );
    }
    
    }
    
    
    // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
    
    
    Foam::functionObjects::runTimePostProcessing::runTimePostProcessing
    
        const Time& runTime,
        const dictionary& dict
    
        fvMeshFunctionObject(name, runTime, dict),
        scene_(runTime, name),
    
    {
        read(dict);
    }
    
    
    // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
    
    
    Foam::functionObjects::runTimePostProcessing::~runTimePostProcessing()
    
    {}
    
    
    // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
    
    
    bool Foam::functionObjects::runTimePostProcessing::read(const dictionary& dict)
    
        Info<< type() << " " << name() << ": reading post-processing data" << endl;
    
    
        scene_.read(dict);
    
        const dictionary& outputDict = dict.subDict("output");
    
        outputDict.readEntry("name", output_.name_);
        outputDict.readEntry("width", output_.width_);
        outputDict.readEntry("height", output_.height_);
    
    
    
        readObjects(dict.subOrEmptyDict("points"), points_);
        readObjects(dict.subOrEmptyDict("lines"), lines_);
        readObjects(dict.subOrEmptyDict("surfaces"), surfaces_);
    
        const dictionary& textDict = dict.subDict("text");
    
    
        for (const entry& dEntry : textDict)
    
                FatalIOErrorInFunction(textDict)
    
                    << "text must be specified in dictionary format"
                    << exit(FatalIOError);
            }
    
    
    bool Foam::functionObjects::runTimePostProcessing::execute()
    
    bool Foam::functionObjects::runTimePostProcessing::write()
    
        Info<< type() << " " << name() <<  " output:" << nl
    
        // (some low-level rendering functionality does not like it)
    
    
        sigFpe::ignore sigFpeHandling; //<- disable in local scope
    
        auto renderWindow = vtkSmartPointer<vtkRenderWindow>::New();
    
        renderWindow->OffScreenRenderingOn();
        renderWindow->SetSize(output_.width_, output_.height_);
    
    
        // Legacy rendering - was deprecated for 8.1.0
        #if (VTK_MAJOR_VERSION < 8) || \
            ((VTK_MAJOR_VERSION == 8) && (VTK_MINOR_VERSION < 2))
        renderWindow->SetAAFrames(10);
    
        renderWindow->SetAlphaBitPlanes(true);
        renderWindow->SetMultiSamples(0);
    //    renderWindow->PolygonSmoothingOn();
    
    
        auto renderer = vtkSmartPointer<vtkRenderer>::New();
    
        scene_.initialise(renderer, output_.name_);
    
        renderWindow->AddRenderer(renderer);
    
        // Add the points
        forAll(points_, i)
        {
            points_[i].addGeometryToScene(0, renderer);
        }
    
        // Add the lines
        forAll(lines_, i)
        {
            lines_[i].addGeometryToScene(0, renderer);
        }
    
        // Add the surfaces
        forAll(surfaces_, i)
        {
            surfaces_[i].addGeometryToScene(0, renderer);
        }
    
    
        // Add the text
        forAll(text_, i)
        {
            text_[i].addGeometryToScene(0, renderer);
        }
    
    
        while (scene_.loop(renderer))
        {
            scalar position = scene_.position();
    
    
                text_[i].updateActors(position);
    
            }
    
            // Update the points
            forAll(points_, i)
            {
                points_[i].updateActors(position);
            }
    
            // Update the lines
            forAll(lines_, i)
            {
                lines_[i].updateActors(position);
            }
    
            // Update the surfaces
            forAll(surfaces_, i)
            {
                surfaces_[i].updateActors(position);
            }
        }
    
        // Clean up
        forAll(text_, i)
        {
            text_[i].clear();
        }
        forAll(points_, i)
        {
            points_[i].clear();
        }
        forAll(lines_, i)
        {
            lines_[i].clear();
        }
        forAll(surfaces_, i)
        {
            surfaces_[i].clear();
        }
    
    
    
        // Instead of relying on the destructor, manually restore the previous
        // SIGFPE state.
        // This is only to avoid compiler complaints about unused variables.
    
        sigFpeHandling.restore();
    
    }
    
    
    // ************************************************************************* //