runTimePostProcessing.C 6.21 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
     \\/     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"
34
#include "sigFpe.H"
35
#include "addToRunTimeSelectionTable.H"
36
37
38
39
40
41
42
43
44
45
46
47

// VTK includes
#include "vtkPolyDataMapper.h"
#include "vtkRenderer.h"
#include "vtkRenderWindow.h"
#include "vtkSmartPointer.h"

#include "vtkLight.h"

// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //

namespace Foam
48
49
{
namespace functionObjects
50
51
{
    defineTypeNameAndDebug(runTimePostProcessing, 0);
52
53
54
55
56
57
58
59

    addToRunTimeSelectionTable
    (
        functionObject,
        runTimePostProcessing,
        dictionary
    );
}
60
61
62
63
64
}


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

65
Foam::functionObjects::runTimePostProcessing::runTimePostProcessing
66
67
(
    const word& name,
68
69
    const Time& runTime,
    const dictionary& dict
70
71
)
:
72
73
    fvMeshFunctionObject(name, runTime, dict),
    scene_(runTime, name),
74
75
76
    points_(),
    lines_(),
    surfaces_(),
77
    text_()
78
79
80
81
82
83
84
{
    read(dict);
}


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

85
Foam::functionObjects::runTimePostProcessing::~runTimePostProcessing()
86
87
88
89
90
{}


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

91
bool Foam::functionObjects::runTimePostProcessing::read(const dictionary& dict)
92
{
93
94
    fvMeshFunctionObject::read(dict);

95
    Info<< type() << " " << name() << ": reading post-processing data" << endl;
96
97
98
99

    scene_.read(dict);

    const dictionary& outputDict = dict.subDict("output");
100
101
102
    outputDict.readEntry("name", output_.name_);
    outputDict.readEntry("width", output_.width_);
    outputDict.readEntry("height", output_.height_);
103
104
105
106
107
108
109


    readObjects(dict.subOrEmptyDict("points"), points_);
    readObjects(dict.subOrEmptyDict("lines"), lines_);
    readObjects(dict.subOrEmptyDict("surfaces"), surfaces_);

    const dictionary& textDict = dict.subDict("text");
110
111

    for (const entry& dEntry : textDict)
112
    {
113
        if (!dEntry.isDict())
114
        {
115
            FatalIOErrorInFunction(textDict)
116
117
118
119
                << "text must be specified in dictionary format"
                << exit(FatalIOError);
        }

120
        text_.append
121
        (
122
123
124
            new runTimePostPro::text
            (
                *this,
125
                dEntry.dict(),
126
127
                scene_.colours()
            )
128
        );
129
130
    }

131
    return true;
132
133
134
}


135
bool Foam::functionObjects::runTimePostProcessing::execute()
136
{
137
    return true;
138
139
140
}


141
bool Foam::functionObjects::runTimePostProcessing::write()
142
143
144
{
    if (!Pstream::master())
    {
145
        return true;
146
147
    }

148
    Info<< type() << " " << name() <<  " output:" << nl
149
150
        << "    Constructing scene" << endl;

151
152

    // Disable any floating point trapping
153
    // (some low-level rendering functionality does not like it)
154
155

    sigFpe::ignore sigFpeHandling; //<- disable in local scope
156

157
    // Initialise render window
158
    auto renderWindow = vtkSmartPointer<vtkRenderWindow>::New();
159
160
    renderWindow->OffScreenRenderingOn();
    renderWindow->SetSize(output_.width_, output_.height_);
161
162
163
164
165

    // Legacy rendering - was deprecated for 8.1.0
    #if (VTK_MAJOR_VERSION < 8) || \
        ((VTK_MAJOR_VERSION == 8) && (VTK_MINOR_VERSION < 2))
    renderWindow->SetAAFrames(10);
166
    #endif
167
168
169
170
    renderWindow->SetAlphaBitPlanes(true);
    renderWindow->SetMultiSamples(0);
//    renderWindow->PolygonSmoothingOn();

171
    auto renderer = vtkSmartPointer<vtkRenderer>::New();
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
    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);
    }

194
195
196
197
198
199
    // Add the text
    forAll(text_, i)
    {
        text_[i].addGeometryToScene(0, renderer);
    }

200
201
202
203
    while (scene_.loop(renderer))
    {
        scalar position = scene_.position();

204
        // Update the text
205
206
        forAll(text_, i)
        {
207
            text_[i].updateActors(position);
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
        }

        // 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);
        }
    }
228

229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
    // 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();
    }

247
248
249
250
251
252

    // Instead of relying on the destructor, manually restore the previous
    // SIGFPE state.
    // This is only to avoid compiler complaints about unused variables.

    sigFpeHandling.restore();
253
254

    return true;
255
256
257
258
}


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