VTK  9.2.6
vtkSurfaceLICInterface.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkSurfaceLICMapper.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
59 #ifndef vtkSurfaceLICInterface_h
60 #define vtkSurfaceLICInterface_h
61 
62 #include "vtkObject.h"
63 #include "vtkRenderingLICOpenGL2Module.h" // For export macro
64 
65 class vtkRenderWindow;
66 class vtkRenderer;
67 class vtkActor;
68 class vtkImageData;
69 class vtkDataObject;
70 class vtkDataArray;
73 class vtkWindow;
74 
75 class VTKRENDERINGLICOPENGL2_EXPORT vtkSurfaceLICInterface : public vtkObject
76 {
77 public:
78  static vtkSurfaceLICInterface* New();
80  void PrintSelf(ostream& os, vtkIndent indent) override;
81 
83 
86  void SetNumberOfSteps(int val);
87  vtkGetMacro(NumberOfSteps, int);
89 
91 
94  void SetStepSize(double val);
95  vtkGetMacro(StepSize, double);
97 
99 
110  void SetNormalizeVectors(int val);
111  vtkBooleanMacro(NormalizeVectors, int);
112  vtkGetMacro(NormalizeVectors, int);
114 
116 
121  void SetMaskOnSurface(int val);
122  vtkBooleanMacro(MaskOnSurface, int);
123  vtkGetMacro(MaskOnSurface, int);
125 
127 
143  void SetMaskThreshold(double val);
144  vtkGetMacro(MaskThreshold, double);
146 
148 
153  void SetMaskColor(double* val);
154  void SetMaskColor(double r, double g, double b)
155  {
156  double rgb[3] = { r, g, b };
157  this->SetMaskColor(rgb);
158  }
159  vtkGetVector3Macro(MaskColor, double);
161 
163 
171  void SetMaskIntensity(double val);
172  vtkGetMacro(MaskIntensity, double);
174 
176 
181  void SetEnhancedLIC(int val);
182  vtkGetMacro(EnhancedLIC, int);
183  vtkBooleanMacro(EnhancedLIC, int);
185 
187 
220  enum
221  {
222  ENHANCE_CONTRAST_OFF = 0,
223  ENHANCE_CONTRAST_LIC = 1,
224  ENHANCE_CONTRAST_COLOR = 3,
225  ENHANCE_CONTRAST_BOTH = 4
226  };
227  void SetEnhanceContrast(int val);
228  vtkGetMacro(EnhanceContrast, int);
230 
232 
248  vtkGetMacro(LowLICContrastEnhancementFactor, double);
249  vtkGetMacro(HighLICContrastEnhancementFactor, double);
250  void SetLowLICContrastEnhancementFactor(double val);
251  void SetHighLICContrastEnhancementFactor(double val);
252  //
253  vtkGetMacro(LowColorContrastEnhancementFactor, double);
254  vtkGetMacro(HighColorContrastEnhancementFactor, double);
255  void SetLowColorContrastEnhancementFactor(double val);
256  void SetHighColorContrastEnhancementFactor(double val);
258 
260 
266  void SetAntiAlias(int val);
267  vtkBooleanMacro(AntiAlias, int);
268  vtkGetMacro(AntiAlias, int);
270 
272 
281  enum
282  {
283  COLOR_MODE_BLEND = 0,
284  COLOR_MODE_MAP
285  };
286  void SetColorMode(int val);
287  vtkGetMacro(ColorMode, int);
289 
291 
300  void SetLICIntensity(double val);
301  vtkGetMacro(LICIntensity, double);
303 
305 
312  void SetMapModeBias(double val);
313  vtkGetMacro(MapModeBias, double);
315 
317 
322  void SetNoiseDataSet(vtkImageData* data);
323  vtkImageData* GetNoiseDataSet();
325 
327 
346  void SetGenerateNoiseTexture(int shouldGenerate);
347  vtkGetMacro(GenerateNoiseTexture, int);
349 
351 
356  enum
357  {
358  NOISE_TYPE_UNIFORM = 0,
359  NOISE_TYPE_GAUSSIAN = 1,
360  NOISE_TYPE_PERLIN = 2
361  };
362  void SetNoiseType(int type);
363  vtkGetMacro(NoiseType, int);
365 
367 
371  void SetNoiseTextureSize(int length);
372  vtkGetMacro(NoiseTextureSize, int);
374 
376 
379  void SetNoiseGrainSize(int val);
380  vtkGetMacro(NoiseGrainSize, int);
382 
384 
390  void SetMinNoiseValue(double val);
391  void SetMaxNoiseValue(double val);
392  vtkGetMacro(MinNoiseValue, double);
393  vtkGetMacro(MaxNoiseValue, double);
395 
397 
401  void SetNumberOfNoiseLevels(int val);
402  vtkGetMacro(NumberOfNoiseLevels, int);
404 
406 
410  void SetImpulseNoiseProbability(double val);
411  vtkGetMacro(ImpulseNoiseProbability, double);
413 
415 
418  void SetImpulseNoiseBackgroundValue(double val);
419  vtkGetMacro(ImpulseNoiseBackgroundValue, double);
421 
423 
426  void SetNoiseGeneratorSeed(int val);
427  vtkGetMacro(NoiseGeneratorSeed, int);
429 
431 
434  enum
435  {
436  COMPOSITE_INPLACE = 0,
437  COMPOSITE_INPLACE_DISJOINT = 1,
438  COMPOSITE_BALANCED = 2,
439  COMPOSITE_AUTO = 3
440  };
441  void SetCompositeStrategy(int val);
442  vtkGetMacro(CompositeStrategy, int);
444 
449  static bool IsSupported(vtkRenderWindow* context);
450 
457  virtual void WriteTimerLog(const char*) {}
458 
462  void ShallowCopy(vtkSurfaceLICInterface* m);
463 
469  virtual void ReleaseGraphicsResources(vtkWindow* win);
470 
474  bool CanRenderSurfaceLIC(vtkActor* actor);
475 
479  void ValidateContext(vtkRenderer* renderer);
480 
487  virtual vtkPainterCommunicator* CreateCommunicator(int);
488 
493  void CreateCommunicator(vtkRenderer*, vtkActor*, vtkDataObject* data);
494 
495  vtkPainterCommunicator* GetCommunicator();
496 
501  void UpdateCommunicator(vtkRenderer* renderer, vtkActor* actor, vtkDataObject* data);
502 
504 
507  void SetHasVectors(bool val);
508  bool GetHasVectors();
510 
514  void InitializeResources();
515 
516  void PrepareForGeometry();
517  void CompletedGeometry();
518  void GatherVectors();
519  void ApplyLIC();
520  void CombineColorsAndLIC();
521  void CopyToScreen();
522 
528  virtual void GetGlobalMinMax(vtkPainterCommunicator*, float&, float&) {}
529 
531 
534  vtkSetMacro(Enable, int);
535  vtkGetMacro(Enable, int);
536  vtkBooleanMacro(Enable, int);
538 
539 protected:
541  ~vtkSurfaceLICInterface() override;
542 
546  void UpdateNoiseImage(vtkRenderWindow* renWin);
547 
549 
552  virtual bool NeedToUpdateCommunicator();
553  bool NeedToRenderGeometry(vtkRenderer* renderer, vtkActor* actor);
554  bool NeedToGatherVectors();
555  bool NeedToComputeLIC();
556  bool NeedToColorLIC();
557  void SetUpdateAll();
559 
560  int Enable;
561 
562  // Unit is a pixel length.
564  double StepSize;
566 
574 
578  double MaskColor[3];
579 
581  double LICIntensity;
582  double MapModeBias;
583 
594 
597 
599 
600 private:
602  void operator=(const vtkSurfaceLICInterface&) = delete;
603 };
604 
605 #endif
void SetMaskColor(double r, double g, double b)
The MaskColor is used on masked fragments.
represents an object (geometry & properties) in a rendered scene
Definition: vtkActor.h:51
abstract base class for most VTK objects
Definition: vtkObject.h:62
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
abstract specification for renderers
Definition: vtkRenderer.h:72
vtkSurfaceLICHelper * Internals
window superclass for vtkRenderWindow
Definition: vtkWindow.h:38
a simple class to control print indentation
Definition: vtkIndent.h:39
topologically and geometrically regular array of data
Definition: vtkImageData.h:53
virtual void WriteTimerLog(const char *)
Methods used for parallel benchmarks.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
virtual void GetGlobalMinMax(vtkPainterCommunicator *, float &, float &)
Get the min/max across all ranks.
create a window for renderers to draw into
public API for surface lic parameters arbitrary geometry.
A communicator that can safely be used inside a painter.
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
general representation of visualization data
Definition: vtkDataObject.h:65
A small collection of noise routines for LIC.