scene.C 11.9 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2015-2016 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 30 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
     \\/     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 "scene.H"
#include "Constant.H"

// VTK includes
#include "vtkCamera.h"
#include "vtkCubeSource.h"
#include "vtkLightKit.h"
#include "vtkPolyDataMapper.h"
#include "vtkPNGWriter.h"
#include "vtkRenderer.h"
#include "vtkRendererCollection.h"
#include "vtkRenderWindow.h"
#include "vtkWindowToImageFilter.h"

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

namespace Foam
{
    template<>
    const char* NamedEnum<scene::modeType, 2>::names[] =
    {
        "static",
        "flightPath"
    };
}

const Foam::NamedEnum<Foam::scene::modeType, 2> modeTypeNames_;


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

void Foam::scene::readCamera(const dictionary& dict)
{
    if (dict.readIfPresent("nFrameTotal", nFrameTotal_))
    {
        if (nFrameTotal_ < 1)
        {
64 65
            FatalIOErrorInFunction(dict)
                << "nFrameTotal must be 1 or greater"
66 67 68 69
                << exit(FatalIOError);
        }
    }

70
    if (dict.readIfPresent("startPosition", startPosition_))
71
    {
72
        if ((startPosition_ < 0) || (startPosition_ > 1))
73
        {
74 75
            FatalIOErrorInFunction(dict)
                << "startPosition must be in the range 0-1"
76 77
                << exit(FatalIOError);
        }
78 79 80 81
        else
        {
            position_ = startPosition_;
        }
82 83 84 85 86 87 88 89 90 91
    }


    dict.lookup("parallelProjection") >> parallelProjection_;

    if (nFrameTotal_ > 1)
    {
        scalar endPosition = dict.lookupOrDefault<scalar>("endPosition", 1);
        if ((endPosition < 0) || (endPosition > 1))
        {
92 93
            FatalIOErrorInFunction(dict)
                << "endPosition must be in the range 0-1"
94 95
                << exit(FatalIOError);
        }
96
        dPosition_ = (endPosition - startPosition_)/scalar(nFrameTotal_ - 1);
97 98 99 100 101 102 103 104 105 106 107 108 109
    }

    mode_ = modeTypeNames_.read(dict.lookup("mode"));

    word coeffsName = modeTypeNames_[mode_] + word("Coeffs");
    const dictionary& coeffs = dict.subDict(coeffsName);

    switch (mode_)
    {
        case mtStatic:
        {
            clipBox_ = boundBox(coeffs.lookup("clipBox"));
            const vector lookDir(vector(coeffs.lookup("lookDir")));
110 111 112 113
            cameraPosition_.reset
            (
                new Function1Types::Constant<point>("position", -lookDir)
            );
114 115 116
            const vector focalPoint(coeffs.lookup("focalPoint"));
            cameraFocalPoint_.reset
            (
117
                new Function1Types::Constant<point>("focalPoint", focalPoint)
118 119
            );
            const vector up(coeffs.lookup("up"));
120
            cameraUp_.reset(new Function1Types::Constant<point>("up", up));
121 122 123 124 125 126
            if (nFrameTotal_ > 1)
            {
                WarningInFunction
                    << "Static mode with nFrames > 1 - please check your setup"
                    << endl;
            }
127 128 129 130 131 132
            break;
        }
        case mtFlightPath:
        {
            cameraPosition_.reset
            (
133
                Function1<vector>::New("position", coeffs).ptr()
134 135 136
            );
            cameraFocalPoint_.reset
            (
137
                Function1<point>::New("focalPoint", coeffs).ptr()
138
            );
139
            cameraUp_.reset(Function1<vector>::New("up", coeffs).ptr());
140 141 142 143
            break;
        }
        default:
        {
144
            FatalErrorInFunction
145 146 147 148 149 150 151
                << "Unhandled enumeration " << modeTypeNames_[mode_]
                << abort(FatalError);
        }
    }

    if (dict.found("zoom"))
    {
152
        cameraZoom_.reset(Function1<scalar>::New("zoom", dict).ptr());
153 154 155
    }
    else
    {
156
        cameraZoom_.reset(new Function1Types::Constant<scalar>("zoom", 1.0));
157 158 159 160
    }

    if (dict.found("viewAngle"))
    {
161
        cameraViewAngle_.reset(Function1<scalar>::New("viewAngle", dict).ptr());
162 163 164
    }
    else
    {
165 166 167 168
        cameraViewAngle_.reset
        (
            new Function1Types::Constant<scalar>("viewAngle", 35.0)
        );
169 170 171 172 173 174 175 176 177 178
    }
}


void Foam::scene::readColours(const dictionary& dict)
{
    const wordList colours = dict.toc();
    forAll(colours, i)
    {
        const word& c = colours[i];
179
        colours_.insert(c, Function1<vector>::New(c, dict).ptr());
180 181 182 183 184 185 186
    }
}


void Foam::scene::initialise(vtkRenderer* renderer, const word& outputName)
{
    currentFrameI_ = 0;
187
    position_ = startPosition_;
188 189 190 191

    outputName_ = outputName;

    // Set the background
192
    const vector backgroundColour = colours_["background"]->value(position_);
193 194 195 196 197 198 199 200 201 202 203
    renderer->SetBackground
    (
        backgroundColour.x(),
        backgroundColour.y(),
        backgroundColour.z()
    );

    // Apply gradient background if "background2" defined
    if (colours_.found("background2"))
    {
        renderer->GradientBackgroundOn();
204
        vector backgroundColour2 = colours_["background2"]->value(position_);
205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224

        renderer->SetBackground2
        (
            backgroundColour2.x(),
            backgroundColour2.y(),
            backgroundColour2.z()
        );
    }

    // Depth peeling
    renderer->SetUseDepthPeeling(true);
    renderer->SetMaximumNumberOfPeels(4);
    renderer->SetOcclusionRatio(0);

    // Set the camera
    vtkSmartPointer<vtkCamera> camera = vtkSmartPointer<vtkCamera>::New();
    camera->SetParallelProjection(parallelProjection_);
    renderer->SetActiveCamera(camera);


225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241
    // Initialise the camera
    const vector up = cameraUp_->value(position_);
    const vector pos = cameraPosition_->value(position_);
    const point focalPoint = cameraFocalPoint_->value(position_);

    camera->SetViewUp(up.x(), up.y(), up.z());
    camera->SetPosition(pos.x(), pos.y(), pos.z());
    camera->SetFocalPoint(focalPoint.x(), focalPoint.y(), focalPoint.z());
    camera->Modified();


    // Add the lights
    vtkSmartPointer<vtkLightKit> lightKit = vtkSmartPointer<vtkLightKit>::New();
    lightKit->AddLightsToRenderer(renderer);


    // For static mode initialise the clip box
242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262
    if (mode_ == mtStatic)
    {
        const point& min = clipBox_.min();
        const point& max = clipBox_.max();
        vtkSmartPointer<vtkCubeSource> clipBox =
            vtkSmartPointer<vtkCubeSource>::New();
        clipBox->SetXLength(max.x() - min.x());
        clipBox->SetYLength(max.y() - min.y());
        clipBox->SetZLength(max.z() - min.z());
        clipBox->SetCenter
        (
            min.x() + 0.5*(max.x() - min.x()),
            min.y() + 0.5*(max.y() - min.y()),
            min.z() + 0.5*(max.z() - min.z())
        );
        vtkSmartPointer<vtkPolyDataMapper> clipMapper =
            vtkSmartPointer<vtkPolyDataMapper>::New();
        clipMapper->SetInputConnection(clipBox->GetOutputPort());

        vtkSmartPointer<vtkActor> clipActor = vtkSmartPointer<vtkActor>::New();
        clipActor->SetMapper(clipMapper);
263
        clipActor->VisibilityOff();
264 265
        renderer->AddActor(clipActor);

266 267
        // Call resetCamera to fit clip box in view
        clipActor->VisibilityOn();
268 269 270 271 272 273
        renderer->ResetCamera();
        clipActor->VisibilityOff();
    }
}


274
void Foam::scene::setCamera(vtkRenderer* renderer) const
275
{
276
    if (mode_ == mtFlightPath)
277
    {
278 279 280
        const vector up = cameraUp_->value(position_);
        const vector pos = cameraPosition_->value(position_);
        const point focalPoint = cameraFocalPoint_->value(position_);
281

282
        vtkCamera* camera = renderer->GetActiveCamera();
283 284 285 286
        camera->SetViewUp(up.x(), up.y(), up.z());
        camera->SetPosition(pos.x(), pos.y(), pos.z());
        camera->SetFocalPoint(focalPoint.x(), focalPoint.y(), focalPoint.z());
        camera->Modified();
287
    }
288

289 290 291 292 293 294
    if (!parallelProjection_)
    {
        // Restore viewAngle (it might be reset by clipping)
        vtkCamera* camera = renderer->GetActiveCamera();
        camera->SetViewAngle(cameraViewAngle_->value(position_));
        camera->Modified();
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
    }
}


Foam::string Foam::scene::frameIndexStr() const
{
    string str = Foam::name(currentFrameI_);
    str.insert(0, 4 - str.length(), '0');

    return str;
}


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

Foam::scene::scene(const objectRegistry& obr, const word& name)
:
    obr_(obr),
    name_(name),
    colours_(),
    mode_(mtStatic),
    cameraPosition_(NULL),
    cameraFocalPoint_(NULL),
    cameraUp_(NULL),
    cameraViewAngle_(NULL),
    clipBox_(),
    parallelProjection_(true),
    nFrameTotal_(1),
323
    startPosition_(0),
324 325 326 327 328 329 330 331 332 333 334 335 336 337 338
    position_(0),
    dPosition_(0),
    currentFrameI_(0),
    outputName_("unknown")
{}


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

Foam::scene::~scene()
{}


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

339
const Foam::HashPtrTable<Foam::Function1<Foam::vector>, Foam::word>&
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
Foam::scene::colours() const
{
    return colours_;
}


Foam::label Foam::scene::frameIndex() const
{
    return currentFrameI_;
}


Foam::scalar Foam::scene::position() const
{
    return position_;
}


void Foam::scene::read(const dictionary& dict)
{
    readCamera(dict.subDict("camera"));
    readColours(dict.subDict("colours"));
}


bool Foam::scene::loop(vtkRenderer* renderer)
{
    static bool initialised = false;
368
    setCamera(renderer);
369 370 371 372 373 374 375

    if (!initialised)
    {
        initialised = true;
        return true;
    }

376 377 378 379
    // Ensure that all objects can be seen without clipping
    // Note: can only be done after all objects have been added!
    renderer->ResetCameraClippingRange();

380 381 382 383 384
    // Save image from last iteration
    saveImage(renderer->GetRenderWindow());

    currentFrameI_++;

385
    position_ = startPosition_ + currentFrameI_*dPosition_;
386

387 388 389 390 391 392
    if (currentFrameI_ < nFrameTotal_)
    {
        return true;
    }
    else
    {
393
        initialised = false;
394 395 396 397 398 399 400 401 402 403 404 405
        return false;
    }
}


void Foam::scene::saveImage(vtkRenderWindow* renderWindow) const
{
    if (!renderWindow)
    {
        return;
    }

406 407 408 409 410 411
    const Time& runTime = obr_.time();

    fileName prefix(Pstream::parRun() ?
        runTime.path()/".."/"postProcessing"/name_/obr_.time().timeName() :
        runTime.path()/"postProcessing"/name_/obr_.time().timeName());

412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441
    mkDir(prefix);

    renderWindow->Render();

    // Set up off-screen rendering
    vtkSmartPointer<vtkWindowToImageFilter> windowToImageFilter =
        vtkSmartPointer<vtkWindowToImageFilter>::New();

    windowToImageFilter->SetInput(renderWindow);

    //// Add alpha channel for transparency
    // windowToImageFilter->SetInputBufferTypeToRGBA();
    windowToImageFilter->SetInputBufferTypeToRGB();

//    windowToImageFilter->ReadFrontBufferOff();
    windowToImageFilter->Update();

    // Save the image
    vtkSmartPointer<vtkPNGWriter> writer = vtkSmartPointer<vtkPNGWriter>::New();
    fileName fName(prefix/outputName_ + '.' + frameIndexStr() + ".png");
    writer->SetFileName(fName.c_str());
    writer->SetInputConnection(windowToImageFilter->GetOutputPort());

    Info<< "    Generating image: " << fName << endl;

    writer->Write();
}


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