fieldVisualisationBase.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 | Copyright (C) 2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 // OpenFOAM includes
27 #include "fieldVisualisationBase.H"
28 #include "runTimePostProcessing.H"
29 
30 // VTK includes
31 #include "vtkArrowSource.h"
32 #include "vtkColorTransferFunction.h"
33 #include "vtkFloatArray.h"
34 #include "vtkGlyph3D.h"
35 #include "vtkLookupTable.h"
36 #include "vtkPointData.h"
37 #include "vtkPolyData.h"
38 #include "vtkPolyDataMapper.h"
39 #include "vtkRenderer.h"
40 #include "vtkScalarBarActor.h"
41 #include "vtkSmartPointer.h"
42 #include "vtkSphereSource.h"
43 #include "vtkTextActor.h"
44 #include "vtkTextProperty.h"
45 
46 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50  template<>
52  {
53  "colour",
54  "field"
55  };
56 
57  template<>
59  {
60  "rainbow",
61  "blueWhiteRed",
62  "fire",
63  "greyscale"
64  };
65 }
66 
69 
72 
73 
74 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
75 
76 void Foam::fieldVisualisationBase::setColourMap(vtkLookupTable* lut) const
77 {
78  label nColours = 256;
79 
80  lut->SetNumberOfColors(nColours);
81 
82  vtkSmartPointer<vtkColorTransferFunction> ctf =
84 
85  switch (colourMap_)
86  {
87  case cmRainbow:
88  {
89  ctf->SetColorSpaceToHSV();
90  ctf->AddRGBPoint(0, 0, 0, 1);
91  ctf->AddRGBPoint(0.5, 0, 1, 0);
92  ctf->AddRGBPoint(1, 1, 0, 0);
93  break;
94  }
95  case cmBlueWhiteRed:
96  {
97  // Values taken from ParaView settings
98  ctf->SetColorSpaceToDiverging();
99  ctf->AddRGBPoint(0.0, 0.231373, 0.298039, 0.752941);
100  ctf->AddRGBPoint(0.5, 0.865003, 0.865003, 0.865003);
101  ctf->AddRGBPoint(1.0, 0.705882, 0.0156863, 0.14902);
102  break;
103  }
104  case cmFire:
105  {
106  // Values taken from ParaView settings
107  ctf->SetColorSpaceToRGB();
108  ctf->AddRGBPoint(0, 0, 0, 0);
109  ctf->AddRGBPoint(0.4, 0.901961, 0, 0);
110  ctf->AddRGBPoint(0.8, 0.901961, 0.901961, 0);
111  ctf->AddRGBPoint(1, 1, 1, 1);
112  break;
113  }
114  case cmGreyscale:
115  {
116  ctf->SetColorSpaceToRGB();
117  ctf->AddRGBPoint(0, 0, 0, 0);
118  ctf->AddRGBPoint(1, 1, 1, 1);
119  break;
120  }
121  }
122 
123 
124  for (label i = 0; i < nColours; i++)
125  {
126  double* c = ctf->GetColor(scalar(i)/scalar(nColours));
127  lut->SetTableValue(i, c[0], c[1], c[2], 1.0);
128  }
129 }
130 
131 
133 (
134  const scalar position,
135  vtkRenderer* renderer,
136  vtkLookupTable* lut
137 ) const
138 {
139  // add scalar bar legend
140  if (!scalarBar_.visible_)
141  {
142  return;
143  }
144 
145  vtkSmartPointer<vtkScalarBarActor> sbar =
147  sbar->SetLookupTable(lut);
148  sbar->SetNumberOfLabels(scalarBar_.numberOfLabels_);
149 
150  const vector textColour = colours_["text"]->value(position);
151 
152  // workaround to supply our own scalarbar title
153  vtkSmartPointer<vtkTextActor> titleActor =
155  sbar->SetTitle(" ");
156  titleActor->SetInput(scalarBar_.title_.c_str());
157  titleActor->GetTextProperty()->SetFontFamilyToArial();
158  titleActor->GetTextProperty()->SetFontSize(3*scalarBar_.fontSize_);
159  titleActor->GetTextProperty()->SetJustificationToCentered();
160  titleActor->GetTextProperty()->SetVerticalJustificationToBottom();
161  titleActor->GetTextProperty()->BoldOn();
162  titleActor->GetTextProperty()->ItalicOff();
163  titleActor->GetTextProperty()->SetColor
164  (
165  textColour[0],
166  textColour[1],
167  textColour[2]
168  );
169  titleActor->GetPositionCoordinate()->
170  SetCoordinateSystemToNormalizedViewport();
171 
172 /*
173  sbar->SetTitle(scalarBar_.title_.c_str());
174  sbar->GetTitleTextProperty()->SetColor
175  (
176  textColour[0],
177  textColour[1],
178  textColour[2]
179  );
180  sbar->GetTitleTextProperty()->SetFontSize(scalarBar_.fontSize_);
181  sbar->GetTitleTextProperty()->ShadowOff();
182  sbar->GetTitleTextProperty()->BoldOn();
183  sbar->GetTitleTextProperty()->ItalicOff();
184 */
185 
186  sbar->GetLabelTextProperty()->SetColor
187  (
188  textColour[0],
189  textColour[1],
190  textColour[2]
191  );
192  sbar->GetLabelTextProperty()->SetFontSize(scalarBar_.fontSize_);
193  sbar->GetLabelTextProperty()->ShadowOff();
194  sbar->GetLabelTextProperty()->BoldOff();
195  sbar->GetLabelTextProperty()->ItalicOff();
196  sbar->SetLabelFormat(scalarBar_.labelFormat_.c_str());
197 
198  sbar->GetPositionCoordinate()->SetCoordinateSystemToNormalizedViewport();
199  sbar->GetPositionCoordinate()->SetValue
200  (
201  scalarBar_.position_.first(),
202  scalarBar_.position_.second()
203  );
204  if (scalarBar_.vertical_)
205  {
206  sbar->SetOrientationToVertical();
207  sbar->SetWidth(0.1);
208  sbar->SetHeight(0.75);
209  sbar->SetTextPositionToSucceedScalarBar();
210  }
211  else
212  {
213  sbar->SetOrientationToHorizontal();
214 
215  // adjustments since not using scalarbar title property
216  sbar->SetWidth(0.75);
217  sbar->SetHeight(0.07);
218  sbar->SetBarRatio(0.5);
219 // sbar->SetHeight(0.1);
220 // sbar->SetTitleRatio(0.01);
221  sbar->SetTextPositionToPrecedeScalarBar();
222  }
223 
224  titleActor->GetPositionCoordinate()->SetValue
225  (
226  scalarBar_.position_.first() + 0.5*sbar->GetWidth(),
227  scalarBar_.position_.second() + sbar->GetHeight()
228  );
229 
230 // sbar->DrawFrameOn();
231 // sbar->DrawBackgroundOn();
232 // sbar->UseOpacityOff();
233 // sbar->VisibilityOff();
234  sbar->VisibilityOn();
235 
236  renderer->AddActor(sbar);
237  renderer->AddActor2D(titleActor);
238 }
239 
240 
242 (
243  const scalar position,
244  const word& colourFieldName,
245  vtkPolyDataMapper* mapper,
246  vtkRenderer* renderer
247 ) const
248 {
249  mapper->InterpolateScalarsBeforeMappingOn();
250 
251  switch (colourBy_)
252  {
253  case cbColour:
254  {
255  mapper->ScalarVisibilityOff();
256  break;
257  }
258  case cbField:
259  {
260  // create look-up table for colours
261  vtkSmartPointer<vtkLookupTable> lut =
263  setColourMap(lut);
264  lut->SetVectorMode(vtkScalarsToColors::MAGNITUDE);
265 
266  // configure the mapper
267  mapper->SelectColorArray(colourFieldName.c_str());
268  mapper->SetScalarRange(range_.first(), range_.second());
269  mapper->SetScalarModeToUsePointFieldData();
270 
271  mapper->SetColorModeToMapScalars();
272  mapper->SetLookupTable(lut);
273  mapper->ScalarVisibilityOn();
274 
275  // add the bar
276  addScalarBar(position, renderer, lut);
277  break;
278  }
279  }
280 
281  mapper->Modified();
282 }
283 
284 
285 
287 (
288  const scalar position,
289  const word& scaleFieldName,
290  const word& colourFieldName,
291  const scalar maxGlyphLength,
292  vtkPolyData* data,
293  vtkActor* actor,
294  vtkRenderer* renderer
295 ) const
296 {
297  vtkSmartPointer<vtkGlyph3D> glyph = vtkSmartPointer<vtkGlyph3D>::New();
298  vtkSmartPointer<vtkPolyDataMapper> glyphMapper =
300  glyphMapper->SetInputConnection(glyph->GetOutputPort());
301 
302  glyph->SetInputData(data);
303  glyph->ScalingOn();
304  bool ok = true;
305 
306  label nComponents =
307  data->GetPointData()->GetArray(scaleFieldName.c_str())
308  ->GetNumberOfComponents();
309 
310  if (nComponents == 1)
311  {
312  vtkSmartPointer<vtkSphereSource> sphere =
314  sphere->SetCenter(0, 0, 0);
315  sphere->SetRadius(0.5);
316 // setting higher resolution slows the rendering significantly
317 // sphere->SetPhiResolution(20);
318 // sphere->SetThetaResolution(20);
319 
320  glyph->SetSourceConnection(sphere->GetOutputPort());
321 
322  if (maxGlyphLength > 0)
323  {
324  double range[2];
325 
326 // can use values to find range
327 // vtkDataArray* values =
328 // data->GetPointData()->GetScalars(scaleFieldName.c_str());
329 // values->GetRange(range);
330 
331  // set range accoding to user-supplied limits
332  range[0] = range_.first();
333  range[1] = range_.second();
334  glyph->ClampingOn();
335  glyph->SetRange(range);
336 
337  // if range[0] != min(value), maxGlyphLength behaviour will not
338  // be correct...
339  glyph->SetScaleFactor(maxGlyphLength);
340  }
341  else
342  {
343  glyph->SetScaleFactor(1);
344  }
345  glyph->SetScaleModeToScaleByScalar();
346  glyph->OrientOff();
347  glyph->SetColorModeToColorByScalar();
348  glyph->SetInputArrayToProcess
349  (
350  0, // scalars
351  0,
352  0,
353  vtkDataObject::FIELD_ASSOCIATION_POINTS,
354  scaleFieldName.c_str()
355  );
356  }
357  else if (nComponents == 3)
358  {
359  vtkSmartPointer<vtkArrowSource> arrow =
361  arrow->SetTipResolution(10);
362  arrow->SetTipRadius(0.1);
363  arrow->SetTipLength(0.35);
364  arrow->SetShaftResolution(10);
365  arrow->SetShaftRadius(0.03);
366 
367  glyph->SetSourceConnection(arrow->GetOutputPort());
368 
369  if (maxGlyphLength > 0)
370  {
371  vtkDataArray* values =
372  data->GetPointData()->GetVectors(scaleFieldName.c_str());
373  double range[6];
374  values->GetRange(range);
375 
376 /*
377  // attempt to set range for vectors...
378  scalar x0 = sqrt(sqr(range_.first())/3.0);
379  scalar x1 = sqrt(sqr(range_.second())/3.0);
380  range[0] = x0;
381  range[1] = x0;
382  range[2] = x0;
383  range[3] = x1;
384  range[4] = x1;
385  range[5] = x1;
386 */
387  glyph->ClampingOn();
388  glyph->SetRange(range);
389  glyph->SetScaleFactor(maxGlyphLength);
390  }
391  else
392  {
393  glyph->SetScaleFactor(1);
394  }
395  glyph->SetScaleModeToScaleByVector();
396  glyph->OrientOn();
397  glyph->SetVectorModeToUseVector();
398  glyph->SetColorModeToColorByVector();
399  glyph->SetInputArrayToProcess
400  (
401  1, // vectors
402  0,
403  0,
404  vtkDataObject::FIELD_ASSOCIATION_POINTS,
405  scaleFieldName.c_str()
406  );
407  }
408  else
409  {
411  << "Glyphs can only be added to " << pTraits<scalar>::typeName
412  << " and " << pTraits<vector>::typeName << " fields. "
413  << " Field " << scaleFieldName << " has " << nComponents
414  << " components" << endl;
415 
416  ok = false;
417  }
418 
419  if (ok)
420  {
421  glyph->Update();
422 
423  setField(position, colourFieldName, glyphMapper, renderer);
424 
425  glyphMapper->Update();
426 
427  actor->SetMapper(glyphMapper);
428 
429  renderer->AddActor(actor);
430  }
431 }
432 
433 
434 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
435 
437 (
438  const runTimePostProcessing& parent,
439  const dictionary& dict,
440  const HashPtrTable<DataEntry<vector>, word>& colours
441 )
442 :
443  parent_(parent),
444  colours_(colours),
445  fieldName_(dict.lookup("fieldName")),
446  colourBy_(cbColour),
447  colourMap_(cmRainbow),
448  range_()
449 {
450  colourBy_ = colourByTypeNames.read(dict.lookup("colourBy"));
451 
452  switch (colourBy_)
453  {
454  case cbColour:
455  {
456  break;
457  }
458  case cbField:
459  {
460  dict.lookup("range") >> range_;
461  break;
462  }
463  }
464 
465  if (dict.found("colourMap"))
466  {
467  colourMap_ = colourMapTypeNames.read(dict.lookup("colourMap"));
468  }
469 
470  const dictionary& sbarDict = dict.subDict("scalarBar");
471  sbarDict.lookup("visible") >> scalarBar_.visible_;
472  if (scalarBar_.visible_)
473  {
474  sbarDict.lookup("vertical") >> scalarBar_.vertical_;
475  sbarDict.lookup("position") >> scalarBar_.position_;
476  sbarDict.lookup("title") >> scalarBar_.title_;
477  sbarDict.lookup("fontSize") >> scalarBar_.fontSize_;
478  sbarDict.lookup("labelFormat") >> scalarBar_.labelFormat_;
479  sbarDict.lookup("numberOfLabels") >> scalarBar_.numberOfLabels_;
480  }
481 }
482 
483 
484 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
485 
487 {}
488 
489 
490 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
491 
494 {
495  return colours_;
496 }
497 
498 
500 {
501  return fieldName_;
502 }
503 
504 
505 // ************************************************************************* //
Foam::fieldVisualisationBase::setColourMap
void setColourMap(vtkLookupTable *lut) const
Set the colour map.
Definition: fieldVisualisationBase.C:76
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::fieldVisualisationBase::setField
void setField(const scalar position, const word &colourFieldName, vtkPolyDataMapper *mapper, vtkRenderer *renderer) const
Set field/configure mapper, add scalar bar.
Definition: fieldVisualisationBase.C:242
Foam::fieldVisualisationBase::cmGreyscale
@ cmGreyscale
Definition: fieldVisualisationBase.H:84
Foam::fieldVisualisationBase::cmRainbow
@ cmRainbow
Definition: fieldVisualisationBase.H:81
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::fieldVisualisationBase::colourMap_
colourMapType colourMap_
Colour map type.
Definition: fieldVisualisationBase.H:132
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::runTimePostProcessing
Function object to generate images during run-time.
Definition: runTimePostProcessing.H:83
Foam::fieldVisualisationBase::cmBlueWhiteRed
@ cmBlueWhiteRed
Definition: fieldVisualisationBase.H:82
Foam::fieldVisualisationBase::~fieldVisualisationBase
virtual ~fieldVisualisationBase()
Destructor.
Definition: fieldVisualisationBase.C:486
Foam::fieldVisualisationBase::addGlyphs
void addGlyphs(const scalar position, const word &scaleFieldName, const word &colourFieldName, const scalar maxGlyphLength, vtkPolyData *data, vtkActor *actor, vtkRenderer *renderer) const
Add glyphs.
Definition: fieldVisualisationBase.C:287
Foam::fieldVisualisationBase::cmFire
@ cmFire
Definition: fieldVisualisationBase.H:83
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::fieldVisualisationBase::colours
const HashPtrTable< DataEntry< vector >, word > & colours() const
Return the colours.
Definition: fieldVisualisationBase.C:493
setField
surfacesMesh setField(triSurfaceToAgglom)
Foam::fieldVisualisationBase::colourMapTypeNames
static const NamedEnum< colourMapType, 4 > colourMapTypeNames
Definition: fieldVisualisationBase.H:87
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
fieldVisualisationBase.H
Foam::fieldVisualisationBase::fieldName
const word & fieldName() const
Return the field name.
Definition: fieldVisualisationBase.C:499
runTimePostProcessing.H
Foam::fieldVisualisationBase::fieldVisualisationBase
fieldVisualisationBase(const fieldVisualisationBase &)
Disallow default bitwise copy construct.
range
scalar range
Definition: LISASMDCalcMethod1.H:12
Foam::HashPtrTable
A HashTable specialization for hashing pointers.
Definition: HashPtrTable.H:50
Foam::Vector< scalar >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::ITstream::read
virtual Istream & read(token &)
Return next token from stream.
Definition: ITstream.C:56
Foam::fieldVisualisationBase::colourByTypeNames
static const NamedEnum< colourByType, 2 > colourByTypeNames
Definition: fieldVisualisationBase.H:77
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::fieldVisualisationBase::addScalarBar
void addScalarBar(const scalar position, vtkRenderer *renderer, vtkLookupTable *lut) const
Add scalar bar to renderer.
Definition: fieldVisualisationBase.C:133
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::DataEntry
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: DataEntry.H:52
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:52
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52