runTimePostProcessing.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2019 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 // OpenFOAM includes
29 #include "runTimePostProcessing.H"
30 #include "dictionary.H"
31 #include "pointData.H"
32 #include "pathline.H"
33 #include "surface.H"
34 #include "text.H"
35 #include "Time.H"
36 #include "sigFpe.H"
37 #include "polySurfaceFields.H"
39 
40 // VTK includes
41 #include "vtkPolyDataMapper.h"
42 #include "vtkRenderer.h"
43 #include "vtkRenderWindow.h"
44 #include "vtkSmartPointer.h"
45 #include "vtkLight.h"
46 
47 #include "vtkDummyController.h"
48 
49 #ifdef FOAM_USING_VTK_MPI
50 # include "vtkMPICommunicator.h"
51 # include "vtkMPIController.h"
52 # include "vtkCompositedSynchronizedRenderers.h"
53 # include "vtkSynchronizedRenderWindows.h"
54 #endif
55 
56 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 namespace functionObjects
61 {
62  defineTypeNameAndDebug(runTimePostProcessing, 0);
63 
65  (
66  functionObject,
67  runTimePostProcessing,
68  dictionary
69  );
70 }
71 }
72 
73 
74 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
75 
76 namespace Foam
77 {
78 
79 template<class Type>
80 static void addGeometryToScene
81 (
82  PtrList<Type>& objects,
83  const scalar position,
84  vtkRenderer* renderer
85 )
86 {
87  for (Type& obj : objects)
88  {
89  if (Pstream::master() || obj.parallel())
90  {
91  obj.addGeometryToScene(position, renderer);
92  }
93  else
94  {
95  obj.addGeometryToScene(position, nullptr);
96  }
97  }
98 }
99 
100 
101 template<class Type>
102 static void updateActors(PtrList<Type>& objects, const scalar position)
103 {
104  for (Type& obj : objects)
105  {
106  obj.updateActors(position);
107  }
108 }
109 
110 
111 template<class Type>
112 static void cleanup(PtrList<Type>& objects)
113 {
114  for (Type& obj : objects)
115  {
116  obj.clear();
117  }
118 }
119 
120 } // End namespace Foam
121 
122 
123 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
124 
125 void Foam::functionObjects::runTimePostProcessing::render
126 (
127  vtkMultiProcessController* controller
128 )
129 {
130  // Some feedback
131  if (controller)
132  {
133  Log << name() << " render (" << controller->GetNumberOfProcesses()
134  << " processes)" << endl;
135  }
136  else
137  {
138  Log << name() << " render" << endl;
139  }
140 
141  // Disable any floating point trapping
142  // (some low-level rendering functionality does not like it)
143 
144  sigFpe::ignore sigFpeHandling; //<- disable in local scope
145 
146 
147  // Normal rendering elements
150 
151  // Initialise render window
152  if (controller || Pstream::master())
153  {
155  renderWindow = vtkSmartPointer<vtkRenderWindow>::New();
156 
157  renderWindow->OffScreenRenderingOn();
158  renderWindow->SetSize(output_.width_, output_.height_);
159 
160  // Legacy rendering - was deprecated for 8.1.0
161  #if (VTK_MAJOR_VERSION < 8) || \
162  ((VTK_MAJOR_VERSION == 8) && (VTK_MINOR_VERSION < 2))
163  renderWindow->SetAAFrames(10);
164  #endif
165  renderWindow->SetAlphaBitPlanes(true);
166  renderWindow->SetMultiSamples(0);
167  // renderWindow->PolygonSmoothingOn();
168 
169  renderWindow->AddRenderer(renderer);
170  }
171 
172 
173  // ---------------------
174  #ifdef FOAM_USING_VTK_MPI
175 
176  // Multi-process synchronization
179 
180  if (controller)
181  {
182  syncRenderers =
184 
185  syncWindows =
187 
188  syncWindows->SetRenderWindow(renderWindow);
189  syncWindows->SetParallelController(controller);
190  syncWindows->SetIdentifier(1);
191 
192  // false = Call Render() manually on each process - don't use RMI
193  syncWindows->SetParallelRendering(true);
194 
195  syncRenderers->SetRenderer(renderer);
196  syncRenderers->SetParallelController(controller);
197  }
198  #endif
199  // ---------------------
200 
201  scene_.initialise(renderer, output_.name_);
202 
203  addGeometryToScene(points_, 0, renderer);
204  addGeometryToScene(lines_, 0, renderer);
205  addGeometryToScene(surfaces_, 0, renderer);
206  addGeometryToScene(text_, 0, renderer);
207 
208  while (scene_.loop(renderer))
209  {
210  const scalar position = scene_.position();
211 
212  updateActors(text_, position);
213  updateActors(points_, position);
214  updateActors(lines_, position);
215  updateActors(surfaces_, position);
216  }
217 
218  // Cleanup
219  cleanup(text_);
220  cleanup(points_);
221  cleanup(lines_);
222  cleanup(surfaces_);
223 
224 
225  // Instead of relying on the destructor, manually restore the previous
226  // SIGFPE state.
227  // This is only to avoid compiler complaints about unused variables.
228 
229  sigFpeHandling.restore();
230 }
231 
232 
233 void Foam::functionObjects::runTimePostProcessing::render
234 (
235  vtkMultiProcessController* controller,
236  void* processData
237 )
238 {
239  reinterpret_cast<runTimePostProcessing*>(processData)->render(controller);
240 }
241 
242 
243 void Foam::functionObjects::runTimePostProcessing::render()
244 {
245  #ifdef FOAM_USING_VTK_MPI
246  if (parallel_)
247  {
248  // Create vtkMPIController if MPI is configured,
249  // vtkThreadedController otherwise.
251  ctrl->Initialize(nullptr, nullptr, 1);
252 
253  ctrl->SetSingleMethod(runTimePostProcessing::render, this);
254  ctrl->SingleMethodExecute();
255 
256  ctrl->Finalize(1);
257  }
258  else
259  #endif
260  {
261  // Normally we would have a fallback controller like this:
262 
263  // if (Pstream::master())
264  // {
265  // auto ctrl = vtkSmartPointer<vtkDummyController>::New();
266  // ctrl->Initialize(nullptr, nullptr, 1);
267  //
268  // ctrl->SetSingleMethod(runTimePostProcessing::render, this);
269  // ctrl->SingleMethodExecute();
270  //
271  // ctrl->Finalize(1);
272  // }
273 
274  // However, this would prevent us from doing any of our own MPI
275  // since this would only be spawned the master.
276 
277  // Instead pass in nullptr for the controller and handling
278  // logic internally.
279 
280  vtkDummyController* dummy = nullptr;
281  render(dummy);
282  }
283 }
284 
285 
286 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
287 
289 (
290  const word& name,
291  const Time& runTime,
292  const dictionary& dict
293 )
294 :
296  output_(),
297  parallel_(false),
298  scene_(runTime, name),
299  points_(),
300  lines_(),
301  surfaces_(),
302  text_()
303 {
304  read(dict);
305 }
306 
307 
308 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
309 
311 {
313 
314  #ifdef FOAM_USING_VTK_MPI
315  parallel_ = (Pstream::parRun() && dict.getOrDefault("parallel", true));
316  #else
317  parallel_ = false;
318  #endif
319 
320  Info<< type() << " " << name() << ": reading post-processing data ("
321  << (parallel_ ? "parallel" : "serial") << " rendering)" << endl;
322 
323  if (dict.getOrDefault("debug", false))
324  {
326  Info<< " debugging on" << endl;
327  }
328 
329  scene_.read(dict);
330 
331  const dictionary& outputDict = dict.subDict("output");
332  outputDict.readEntry("name", output_.name_);
333  outputDict.readEntry("width", output_.width_);
334  outputDict.readEntry("height", output_.height_);
335 
336  readObjects(dict.subOrEmptyDict("points"), points_);
337  readObjects(dict.subOrEmptyDict("lines"), lines_);
338  readObjects(dict.subOrEmptyDict("surfaces"), surfaces_);
339 
340  const dictionary& textDict = dict.subDict("text");
341 
342  for (const entry& dEntry : textDict)
343  {
344  if (!dEntry.isDict())
345  {
346  FatalIOErrorInFunction(textDict)
347  << textDict.dictName()
348  << "text must be specified in dictionary format"
349  << exit(FatalIOError);
350  }
351 
352  const dictionary& objectDict = dEntry.dict();
353 
354  text_.append
355  (
357  (
358  *this,
359  objectDict,
360  scene_.colours()
361  )
362  );
363  }
364 
365  return true;
366 }
367 
368 
370 {
371  return true;
372 }
373 
374 
376 {
377  render();
378  return true;
379 }
380 
381 
382 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
surface.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
polySurfaceFields.H
Fields for polySurface.
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:414
pathline.H
Foam::addGeometryToScene
static void addGeometryToScene(PtrList< Type > &objects, const scalar position, vtkRenderer *renderer)
Definition: runTimePostProcessing.C:81
vtkSmartPointer
Definition: runTimePostProcessing.H:148
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::functionObjects::runTimePostPro::scene::loop
bool loop(vtkRenderer *renderer)
Main control loop.
Definition: scene.C:350
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
Foam::functionObjects::runTimePostProcessing::execute
virtual bool execute()
Execute, currently does nothing.
Definition: runTimePostProcessing.C:369
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, add, dictionary)
Foam::updateActors
static void updateActors(PtrList< Type > &objects, const scalar position)
Definition: runTimePostProcessing.C:102
Foam::cleanup
static void cleanup(PtrList< Type > &objects)
Definition: runTimePostProcessing.C:112
Foam::functionObjects::runTimePostPro::scene::initialise
void initialise(vtkRenderer *renderer, const word &outputName)
Definition: scene.C:148
Foam::functionObjects::runTimePostProcessing::runTimePostProcessing
runTimePostProcessing(const word &name, const Time &runTime, const dictionary &dict)
Construct from dictionary.
Definition: runTimePostProcessing.C:289
text.H
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:314
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
Foam::functionObjects::runTimePostPro::text
Define text element for runTimePostProcessing.
Definition: text.H:175
Foam::functionObjects::runTimePostPro::geometryBase::debug
static int debug
Debug switch.
Definition: geometryBase.H:160
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::functionObjects::runTimePostProcessing::read
virtual bool read(const dictionary &dict)
Read the post-processing controls.
Definition: runTimePostProcessing.C:310
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:166
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
Time.H
runTimePostProcessing.H
pointData.H
Foam::functionObject::name
const word & name() const
Return the name of this functionObject.
Definition: functionObject.C:131
Foam::functionObjects::runTimePostProcessing::write
virtual bool write()
Write.
Definition: runTimePostProcessing.C:375
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
dictionary.H
Foam::functionObjects::runTimePostPro::scene::position
scalar position() const
Return the current position (in range 0-1)
Definition: scene.C:334
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
Log
#define Log
Report write to Foam::Info if the local log switch is true.
Definition: messageStream.H:332