Newer
Older
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2015-2016 OpenCFD Ltd.
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
-------------------------------------------------------------------------------
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)
{
FatalIOErrorInFunction(dict)
<< "nFrameTotal must be 1 or greater"
<< exit(FatalIOError);
}
}
if (dict.readIfPresent("startPosition", startPosition_))
if ((startPosition_ < 0) || (startPosition_ > 1))
FatalIOErrorInFunction(dict)
<< "startPosition must be in the range 0-1"
<< exit(FatalIOError);
}
else
{
position_ = startPosition_;
}
}
dict.lookup("parallelProjection") >> parallelProjection_;
if (nFrameTotal_ > 1)
{
scalar endPosition = dict.lookupOrDefault<scalar>("endPosition", 1);
if ((endPosition < 0) || (endPosition > 1))
{
FatalIOErrorInFunction(dict)
<< "endPosition must be in the range 0-1"
<< exit(FatalIOError);
}
dPosition_ = (endPosition - startPosition_)/scalar(nFrameTotal_ - 1);
}
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")));
Andrew Heather
committed
cameraPosition_.reset
(
new Function1Types::Constant<point>("position", -lookDir)
);
const vector focalPoint(coeffs.lookup("focalPoint"));
cameraFocalPoint_.reset
(
Andrew Heather
committed
new Function1Types::Constant<point>("focalPoint", focalPoint)
);
const vector up(coeffs.lookup("up"));
Andrew Heather
committed
cameraUp_.reset(new Function1Types::Constant<point>("up", up));
break;
}
case mtFlightPath:
{
cameraPosition_.reset
(
Andrew Heather
committed
Function1<vector>::New("position", coeffs).ptr()
);
cameraFocalPoint_.reset
(
Andrew Heather
committed
Function1<point>::New("focalPoint", coeffs).ptr()
Andrew Heather
committed
cameraUp_.reset(Function1<vector>::New("up", coeffs).ptr());
break;
}
default:
{
<< "Unhandled enumeration " << modeTypeNames_[mode_]
<< abort(FatalError);
}
}
if (dict.found("viewAngle"))
{
Andrew Heather
committed
cameraViewAngle_.reset(Function1<scalar>::New("viewAngle", dict).ptr());
}
else
{
Andrew Heather
committed
cameraViewAngle_.reset
(
new Function1Types::Constant<scalar>("viewAngle", 35.0)
);
}
}
void Foam::scene::readColours(const dictionary& dict)
{
const wordList colours = dict.toc();
forAll(colours, i)
{
const word& c = colours[i];
Andrew Heather
committed
colours_.insert(c, Function1<vector>::New(c, dict).ptr());
}
}
void Foam::scene::initialise(vtkRenderer* renderer, const word& outputName)
{
currentFrameI_ = 0;
position_ = startPosition_;
outputName_ = outputName;
// Set the background
const vector backgroundColour = colours_["background"]->value(position_);
renderer->SetBackground
(
backgroundColour.x(),
backgroundColour.y(),
backgroundColour.z()
);
// Apply gradient background if "background2" defined
if (colours_.found("background2"))
{
renderer->GradientBackgroundOn();
vector backgroundColour2 = colours_["background2"]->value(position_);
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);
// 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
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);
clipActor->VisibilityOff();
renderer->AddActor(clipActor);
// Call resetCamera to fit clip box in view
clipActor->VisibilityOn();
renderer->ResetCamera();
clipActor->VisibilityOff();
}
}
void Foam::scene::setCamera(vtkRenderer* renderer) const
if (mode_ == mtFlightPath)
const vector up = cameraUp_->value(position_);
const vector pos = cameraPosition_->value(position_);
const point focalPoint = cameraFocalPoint_->value(position_);
vtkCamera* camera = renderer->GetActiveCamera();
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();
}
if (!parallelProjection_)
{
// Restore viewAngle (it might be reset by clipping)
vtkCamera* camera = renderer->GetActiveCamera();
camera->SetViewAngle(cameraViewAngle_->value(position_));
camera->Modified();
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
}
}
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),
position_(0),
dPosition_(0),
currentFrameI_(0),
outputName_("unknown")
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::scene::~scene()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Andrew Heather
committed
const Foam::HashPtrTable<Foam::Function1<Foam::vector>, Foam::word>&
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
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;
setCamera(renderer);
if (!initialised)
{
initialised = true;
return true;
}
// Ensure that all objects can be seen without clipping
// Note: can only be done after all objects have been added!
renderer->ResetCameraClippingRange();
// Save image from last iteration
saveImage(renderer->GetRenderWindow());
currentFrameI_++;
position_ = startPosition_ + currentFrameI_*dPosition_;
if (currentFrameI_ < nFrameTotal_)
{
return true;
}
else
{
initialised = false;
return false;
}
}
void Foam::scene::saveImage(vtkRenderWindow* renderWindow) const
{
if (!renderWindow)
{
return;
}
const Time& runTime = obr_.time();
fileName prefix(Pstream::parRun() ?
runTime.path()/".."/"postProcessing"/name_/obr_.time().timeName() :
runTime.path()/"postProcessing"/name_/obr_.time().timeName());
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
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();
}
// ************************************************************************* //