VTK  9.3.0
vtkHyperTreeGridMapper.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
23#ifndef vtkHyperTreeGridMapper_h
24#define vtkHyperTreeGridMapper_h
25
26#include "vtkMapper.h"
27#include "vtkSetGet.h" // Get macro
28#include "vtkSmartPointer.h" // For vtkSmartPointer
29
30#include "vtkRenderingHyperTreeGridModule.h" // For export macro
31
32VTK_ABI_NAMESPACE_BEGIN
35class vtkPolyData;
37class vtkRenderWindow;
38class vtkRenderer;
39
40class VTKRENDERINGHYPERTREEGRID_EXPORT vtkHyperTreeGridMapper : public vtkMapper
41{
42public:
45 void PrintSelf(ostream& os, vtkIndent indent) override;
46
48
59 using Superclass::SetInputConnection;
60 void SetInputDataObject(int port, vtkDataObject* input) override;
61 void SetInputDataObject(vtkDataObject* input) override;
63
65
70 double* GetBounds() override;
71 void GetBounds(double bounds[6]) override;
73
75
80 vtkGetMacro(UseAdaptiveDecimation, bool);
81 vtkSetMacro(UseAdaptiveDecimation, bool);
82 vtkBooleanMacro(UseAdaptiveDecimation, bool);
84
91 void Render(vtkRenderer* ren, vtkActor* act) override;
92
98 int FillInputPortInformation(int port, vtkInformation* info) override;
99
100protected:
103
109
110 // In 2D mode, these variables control the mapper oprimisations
111 bool UseAdaptiveDecimation = false;
112
113 // render the extracted surface,
114 // need to be created in device specific subclass
116
117 // Internal object to render
119
120private:
122 void operator=(const vtkHyperTreeGridMapper&) = delete;
123};
124
125VTK_ABI_NAMESPACE_END
126#endif
represents an object (geometry & properties) in a rendered scene
Definition vtkActor.h:41
abstract superclass for composite (multi-block or AMR) datasets
general representation of visualization data
map vtkHyperTreeGrid to graphics primitives
void SetInputDataObject(int port, vtkDataObject *input) override
Sets the data-object as an input on the given port index.
void SetInputDataObject(vtkDataObject *input) override
double * GetBounds() override
For this mapper, the bounds correspond to the output for the internal surface filter which may be res...
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
vtkSmartPointer< vtkPolyDataMapper > Mapper
void GetBounds(double bounds[6]) override
For this mapper, the bounds correspond to the output for the internal surface filter which may be res...
void Render(vtkRenderer *ren, vtkActor *act) override
Use the internal PolyData Mapper to do the rendering of the HTG transformed by the current SurfaceFil...
~vtkHyperTreeGridMapper() override
vtkSmartPointer< vtkCompositeDataSet > UpdateWithDecimation(vtkCompositeDataSet *htg, vtkRenderer *ren)
Generate a new composite were each leave is decimated if required.
static vtkHyperTreeGridMapper * New()
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkSmartPointer< vtkCompositeDataSet > Input
A dataset containing a grid of vtkHyperTree instances arranged as a rectilinear grid.
a simple class to control print indentation
Definition vtkIndent.h:29
Store vtkAlgorithm input/output information.
abstract class specifies interface to map data to graphics primitives
Definition vtkMapper.h:77
map vtkPolyData to graphics primitives
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:80
create a window for renderers to draw into
abstract specification for renderers
Definition vtkRenderer.h:59
Hold a reference to a vtkObjectBase instance.