VTK  9.6.20260201
vtkMapper.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
59
60#ifndef vtkMapper_h
61#define vtkMapper_h
62
63#include "vtkAbstractMapper3D.h"
64#include "vtkRenderingCoreModule.h" // For export macro
65#include "vtkSmartPointer.h" // needed for vtkSmartPointer.
66#include "vtkSystemIncludes.h" // For VTK_COLOR_MODE_DEFAULT and _MAP_SCALARS
67#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
68#include <vector> // for method args
69
70#define VTK_RESOLVE_OFF 0
71#define VTK_RESOLVE_POLYGON_OFFSET 1
72#define VTK_RESOLVE_SHIFT_ZBUFFER 2
73
74#define VTK_GET_ARRAY_BY_ID 0
75#define VTK_GET_ARRAY_BY_NAME 1
76
77#define VTK_MATERIALMODE_DEFAULT 0
78#define VTK_MATERIALMODE_AMBIENT 1
79#define VTK_MATERIALMODE_DIFFUSE 2
80#define VTK_MATERIALMODE_AMBIENT_AND_DIFFUSE 3
81
82VTK_ABI_NAMESPACE_BEGIN
83class vtkActor;
84class vtkDataSet;
85class vtkDataObject;
86class vtkFloatArray;
88class vtkImageData;
89class vtkProp;
90class vtkRenderer;
92class vtkSelection;
94class vtkWindow;
95
96class VTKRENDERINGCORE_EXPORT VTK_MARSHALAUTO vtkMapper : public vtkAbstractMapper3D
97{
98public:
100 void PrintSelf(ostream& os, vtkIndent indent) override;
101
106
112
117 virtual void Render(vtkRenderer* ren, vtkActor* a) = 0;
118
125
127
133
139
141
146 vtkBooleanMacro(ScalarVisibility, vtkTypeBool);
148
150
156 vtkSetMacro(Static, vtkTypeBool);
157 vtkGetMacro(Static, vtkTypeBool);
158 vtkBooleanMacro(Static, vtkTypeBool);
160
162
174 vtkSetMacro(ColorMode, int);
175 vtkGetMacro(ColorMode, int);
180
184 const char* GetColorModeAsString();
185
187
197
199
211
213
218 vtkSetVector2Macro(ScalarRange, double);
219 vtkGetVectorMacro(ScalarRange, double, 2);
221
234
235 // When ScalarMode is set to use Field Data (ScalarModeToFieldData),
236 // you must call SelectColorArray to choose the field data array to
237 // be used to color cells. In this mode, the default behavior is to
238 // treat the field data tuples as being associated with cells. If
239 // the poly data contains triangle strips, the array is expected to
240 // contain the cell data for each mini-cell formed by any triangle
241 // strips in the poly data as opposed to treating them as a single
242 // tuple that applies to the entire strip. This mode can also be
243 // used to color the entire poly data by a single color obtained by
244 // mapping the tuple at a given index in the field data array
245 // through the color map. Use SetFieldDataTupleId() to specify
246 // the tuple index.
247 vtkSetMacro(ScalarMode, int);
248 vtkGetMacro(ScalarMode, int);
261
263
268 void SelectColorArray(int arrayNum);
269 void SelectColorArray(const char* arrayName);
271
272 // When ScalarMode is set to UseFieldData, set the index of the
273 // tuple by which to color the entire data set. By default, the
274 // index is -1, which means to treat the field data array selected
275 // with SelectColorArray as having a scalar value for each cell.
276 // Indices of 0 or higher mean to use the tuple at the given index
277 // for coloring the entire data set.
280
282
287 void ColorByArrayComponent(int arrayNum, int component);
288 void ColorByArrayComponent(const char* arrayName, int component);
290
294 vtkGetStringMacro(ArrayName);
295 vtkSetStringMacro(ArrayName);
296 vtkGetMacro(ArrayId, int);
297 vtkSetMacro(ArrayId, int);
298 vtkGetMacro(ArrayAccessMode, int);
299 vtkSetMacro(ArrayAccessMode, int);
300 vtkGetMacro(ArrayComponent, int);
301 vtkSetMacro(ArrayComponent, int);
302
307
309
319 static void SetResolveCoincidentTopology(int val);
331
332
334
339 static void SetResolveCoincidentTopologyPolygonOffsetParameters(double factor, double units);
340 static void GetResolveCoincidentTopologyPolygonOffsetParameters(double& factor, double& units);
342
344
349 void GetRelativeCoincidentTopologyPolygonOffsetParameters(double& factor, double& units);
351
353
358 static void SetResolveCoincidentTopologyLineOffsetParameters(double factor, double units);
359 static void GetResolveCoincidentTopologyLineOffsetParameters(double& factor, double& units);
361
363
367 void SetRelativeCoincidentTopologyLineOffsetParameters(double factor, double units);
368 void GetRelativeCoincidentTopologyLineOffsetParameters(double& factor, double& units);
370
372
380
382
389
391
395 void GetCoincidentTopologyPolygonOffsetParameters(double& factor, double& units);
396 void GetCoincidentTopologyLineOffsetParameters(double& factor, double& units);
399
401
411
413
417 static void SetResolveCoincidentTopologyZShift(double val);
420
425 double* GetBounds() VTK_SIZEHINT(6) override;
426 void GetBounds(double bounds[6]) override { this->vtkAbstractMapper3D::GetBounds(bounds); }
427
433 void SetRenderTime(double time) { this->RenderTime = time; }
434 vtkGetMacro(RenderTime, double);
435
441
443
449 vtkDataSet* GetDataSetInput() { return this->GetInput(); }
451
453
460 virtual vtkUnsignedCharArray* MapScalars(double alpha);
461 virtual vtkUnsignedCharArray* MapScalars(double alpha, int& cellFlag);
462 virtual vtkUnsignedCharArray* MapScalars(vtkDataSet* input, double alpha);
463 virtual vtkUnsignedCharArray* MapScalars(vtkDataSet* input, double alpha, int& cellFlag);
465
467
472 virtual bool HasOpaqueGeometry();
475
482 virtual bool GetSupportsSelection() { return false; }
483
489 std::vector<unsigned int>& /* pixeloffsets */, vtkProp* /* prop */)
490 {
491 }
492
502
508
513
518
523
525
529 vtkGetObjectMacro(Selection, vtkSelection);
532
537 vtkScalarsToColors* lkup, int colorMode);
538
539protected:
541 ~vtkMapper() override;
542
543 // color mapped colors
545
546 // Use texture coordinates for coloring.
548 // Coordinate for each point.
550 // 1D ColorMap used for the texture image.
552 void MapScalarsToTexture(vtkAbstractArray* scalars, double alpha);
553
557 double ScalarRange[2];
559
562
564
565 // for coloring by a component of a field data array
570
571 // If coloring by field data, which tuple to use to color the entire
572 // data set. If -1, treat array values as cell data.
574
576
582
584
585private:
586 vtkMapper(const vtkMapper&) = delete;
587 void operator=(const vtkMapper&) = delete;
588};
589
590VTK_ABI_NAMESPACE_END
591#endif
Abstract superclass for all arrays.
virtual double * GetBounds()=0
Return bounding box (array of six doubles) of data expressed as (xmin,xmax, ymin,ymax,...
represents an object (geometry & properties) in a rendered scene
Definition vtkActor.h:151
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
dynamic, self-adjusting array of float
topologically and geometrically regular array of data
a simple class to control print indentation
Definition vtkIndent.h:108
virtual bool HasOpaqueGeometry()
Some introspection on the type of data the mapper will render used by props to determine if they shou...
const char * GetScalarModeAsString()
Return the method for obtaining scalar data.
~vtkMapper() override
virtual int CanUseTextureMapForColoring(vtkDataObject *input)
Returns if we can use texture maps for scalar coloring.
void GetCoincidentTopologyLineOffsetParameters(double &factor, double &units)
Get the net parameters for handling coincident topology obtained by summing the global values with th...
static int GetResolveCoincidentTopology()
Set/Get a global flag that controls whether coincident topology (e.g., a line on top of a polygon) is...
void SelectColorArray(const char *arrayName)
When ScalarMode is set to UsePointFieldData or UseCellFieldData, you can specify which array to use f...
void SetScalarModeToUsePointFieldData()
Definition vtkMapper.h:252
vtkScalarsToColors * LookupTable
Definition vtkMapper.h:554
double CoincidentPolygonFactor
Definition vtkMapper.h:577
vtkDataSet * GetDataSetInput()
Get the input to this mapper as a vtkDataSet, instead of as a more specialized data type that the sub...
Definition vtkMapper.h:449
void SetColorModeToMapScalars()
default (ColorModeToDefault), unsigned char scalars are treated as colors, and NOT mapped through the...
Definition vtkMapper.h:177
int ColorMode
Definition vtkMapper.h:560
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void SetScalarModeToUseCellData()
Definition vtkMapper.h:251
void SetLookupTable(vtkScalarsToColors *lut)
Specify a lookup table for the mapper to use.
virtual vtkUnsignedCharArray * MapScalars(double alpha)
Map the scalars (if there are any scalars and ScalarVisibility is on) through the lookup table,...
vtkTypeBool InterpolateScalarsBeforeMapping
Definition vtkMapper.h:547
void SetColorModeToDefault()
default (ColorModeToDefault), unsigned char scalars are treated as colors, and NOT mapped through the...
Definition vtkMapper.h:176
void SelectColorArray(int arrayNum)
When ScalarMode is set to UsePointFieldData or UseCellFieldData, you can specify which array to use f...
vtkTypeBool UseLookupTableScalarRange
Definition vtkMapper.h:558
void SetScalarModeToDefault()
Definition vtkMapper.h:249
int ScalarMode
Definition vtkMapper.h:561
vtkSelection * Selection
Definition vtkMapper.h:583
virtual bool GetSupportsSelection()
WARNING: INTERNAL METHOD - NOT INTENDED FOR GENERAL USE DO NOT USE THIS METHOD OUTSIDE OF THE RENDERI...
Definition vtkMapper.h:482
static void GetResolveCoincidentTopologyLineOffsetParameters(double &factor, double &units)
Used to set the line offset scale factor and units.
void GetCoincidentTopologyPolygonOffsetParameters(double &factor, double &units)
Get the net parameters for handling coincident topology obtained by summing the global values with th...
virtual void CreateDefaultLookupTable()
Create default lookup table.
double RenderTime
Definition vtkMapper.h:563
void SetScalarModeToUseCellFieldData()
Definition vtkMapper.h:256
void MapScalarsToTexture(vtkAbstractArray *scalars, double alpha)
vtkImageData * GetColorTextureMap()
Provide read access to the color texture array.
static void GetResolveCoincidentTopologyPointOffsetParameter(double &units)
Used to set the point offset value Used when ResolveCoincidentTopology is set to PolygonOffset.
static void SetResolveCoincidentTopologyToDefault()
Set/Get a global flag that controls whether coincident topology (e.g., a line on top of a polygon) is...
virtual bool HasTranslucentPolygonalGeometry()
Some introspection on the type of data the mapper will render used by props to determine if they shou...
vtkImageData * ColorTextureMap
Definition vtkMapper.h:551
vtkTimeStamp BuildTime
Definition vtkMapper.h:556
void ColorByArrayComponent(int arrayNum, int component)
Legacy: These methods used to be used to specify the array component.
void GetRelativeCoincidentTopologyPointOffsetParameter(double &units)
Used to set the point offset value relative to the global Used when ResolveCoincidentTopology is set ...
void SetColorModeToDirectScalars()
default (ColorModeToDefault), unsigned char scalars are treated as colors, and NOT mapped through the...
Definition vtkMapper.h:178
vtkTypeBool Static
Definition vtkMapper.h:575
virtual vtkUnsignedCharArray * MapScalars(double alpha, int &cellFlag)
Map the scalars (if there are any scalars and ScalarVisibility is on) through the lookup table,...
void SetScalarModeToUseFieldData()
Definition vtkMapper.h:260
static double GetResolveCoincidentTopologyZShift()
Used to set the z-shift if ResolveCoincidentTopology is set to ShiftZBuffer.
virtual void SetColorMode(int)
default (ColorModeToDefault), unsigned char scalars are treated as colors, and NOT mapped through the...
static void SetResolveCoincidentTopologyPolygonOffsetParameters(double factor, double units)
Used to set the polygon offset scale factor and units.
int ArrayAccessMode
Definition vtkMapper.h:569
void SetRelativeCoincidentTopologyLineOffsetParameters(double factor, double units)
Used to set the line offset values relative to the global Used when ResolveCoincidentTopology is set ...
int ArrayId
Definition vtkMapper.h:566
int ArrayComponent
Definition vtkMapper.h:568
void SetRenderTime(double time)
This instance variable is used by vtkLODActor to determine which mapper to use.
Definition vtkMapper.h:433
static void SetResolveCoincidentTopology(int val)
Set/Get a global flag that controls whether coincident topology (e.g., a line on top of a polygon) is...
void GetRelativeCoincidentTopologyLineOffsetParameters(double &factor, double &units)
Used to set the line offset values relative to the global Used when ResolveCoincidentTopology is set ...
vtkUnsignedCharArray * Colors
Definition vtkMapper.h:544
static void SetResolveCoincidentTopologyToOff()
Set/Get a global flag that controls whether coincident topology (e.g., a line on top of a polygon) is...
Definition vtkMapper.h:322
virtual void Render(vtkRenderer *ren, vtkActor *a)=0
Method initiates the mapping process.
vtkTypeBool ScalarVisibility
Definition vtkMapper.h:555
vtkMTimeType GetMTime() override
Overload standard modified time function.
vtkScalarsToColors * GetLookupTable()
Specify a lookup table for the mapper to use.
virtual vtkUnsignedCharArray * MapScalars(vtkDataSet *input, double alpha, int &cellFlag)
Map the scalars (if there are any scalars and ScalarVisibility is on) through the lookup table,...
vtkIdType FieldDataTupleId
Definition vtkMapper.h:573
static vtkSmartPointer< vtkImageData > BuildColorTextureImage(vtkScalarsToColors *lkup, int colorMode)
Create an image of the lookup table lkup.
double CoincidentPointOffset
Definition vtkMapper.h:581
static void SetResolveCoincidentTopologyToPolygonOffset()
Set/Get a global flag that controls whether coincident topology (e.g., a line on top of a polygon) is...
Definition vtkMapper.h:323
static void SetResolveCoincidentTopologyToShiftZBuffer()
Set/Get a global flag that controls whether coincident topology (e.g., a line on top of a polygon) is...
Definition vtkMapper.h:327
void ColorByArrayComponent(const char *arrayName, int component)
Legacy: These methods used to be used to specify the array component.
void GetRelativeCoincidentTopologyPolygonOffsetParameters(double &factor, double &units)
Used to set the polygon offset values relative to the global Used when ResolveCoincidentTopology is s...
static int GetResolveCoincidentTopologyPolygonOffsetFaces()
Used when ResolveCoincidentTopology is set to PolygonOffset.
vtkFloatArray * GetColorCoordinates()
Provide read access to the color texture coordinate array.
double * GetBounds() override
Return bounding box (array of six doubles) of data expressed as (xmin,xmax, ymin,ymax,...
virtual void ProcessSelectorPixelBuffers(vtkHardwareSelector *, std::vector< unsigned int > &, vtkProp *)
allows a mapper to update a selections color buffers Called from a prop which in turn is called from ...
Definition vtkMapper.h:488
double CoincidentLineOffset
Definition vtkMapper.h:580
void SetRelativeCoincidentTopologyPointOffsetParameter(double units)
Used to set the point offset value relative to the global Used when ResolveCoincidentTopology is set ...
static void SetResolveCoincidentTopologyPointOffsetParameter(double units)
Used to set the point offset value Used when ResolveCoincidentTopology is set to PolygonOffset.
vtkFloatArray * ColorCoordinates
Definition vtkMapper.h:549
void SetRelativeCoincidentTopologyPolygonOffsetParameters(double factor, double units)
Used to set the polygon offset values relative to the global Used when ResolveCoincidentTopology is s...
void ShallowCopy(vtkAbstractMapper *m) override
Make a shallow copy of this mapper.
void SetScalarModeToUsePointData()
Definition vtkMapper.h:250
vtkUnsignedCharArray * GetColorMapColors()
Provide read access to the color array.
static void SetResolveCoincidentTopologyLineOffsetParameters(double factor, double units)
Used to set the line offset scale factor and units.
vtkDataSet * GetInput()
Get the input as a vtkDataSet.
virtual void SetSelection(vtkSelection *)
Set/Get selection used to display particular points or cells in a second pass.
static void SetResolveCoincidentTopologyPolygonOffsetFaces(int faces)
Used when ResolveCoincidentTopology is set to PolygonOffset.
static void GetResolveCoincidentTopologyPolygonOffsetParameters(double &factor, double &units)
Used to set the polygon offset scale factor and units.
void GetCoincidentTopologyPointOffsetParameter(double &units)
Get the net parameters for handling coincident topology obtained by summing the global values with th...
double CoincidentLineFactor
Definition vtkMapper.h:579
const char * GetColorModeAsString()
Return the method of coloring scalar data.
static void SetResolveCoincidentTopologyZShift(double val)
Used to set the z-shift if ResolveCoincidentTopology is set to ShiftZBuffer.
virtual vtkUnsignedCharArray * MapScalars(vtkDataSet *input, double alpha)
Map the scalars (if there are any scalars and ScalarVisibility is on) through the lookup table,...
double CoincidentPolygonOffset
Definition vtkMapper.h:578
virtual void SetScalarMode(int)
Control how the filter works with scalar point data and cell attribute data.
void ClearColorArrays()
Call to force a rebuild of color result arrays on next MapScalars.
char * ArrayName
Definition vtkMapper.h:567
double ScalarRange[2]
Definition vtkMapper.h:557
void ReleaseGraphicsResources(vtkWindow *) override
Release any graphics resources that are being consumed by this mapper.
Definition vtkMapper.h:124
abstract superclass for all actors, volumes and annotations
Definition vtkProp.h:69
abstract specification for renderers
Superclass for mapping scalar values to colors.
data object that represents a "selection" in VTK.
Hold a reference to a vtkObjectBase instance.
record modification and/or execution time
dynamic, self-adjusting array of unsigned char
window superclass for vtkRenderWindow
Definition vtkWindow.h:61
int vtkTypeBool
Definition vtkABI.h:64
#define VTK_SCALAR_MODE_DEFAULT
#define VTK_SCALAR_MODE_USE_POINT_DATA
#define VTK_SCALAR_MODE_USE_FIELD_DATA
#define VTK_SCALAR_MODE_USE_CELL_DATA
#define VTK_SCALAR_MODE_USE_CELL_FIELD_DATA
#define VTK_SCALAR_MODE_USE_POINT_FIELD_DATA
#define VTK_RESOLVE_POLYGON_OFFSET
Definition vtkMapper.h:71
#define VTK_RESOLVE_OFF
Definition vtkMapper.h:70
#define VTK_RESOLVE_SHIFT_ZBUFFER
Definition vtkMapper.h:72
#define VTK_COLOR_MODE_MAP_SCALARS
#define VTK_COLOR_MODE_DEFAULT
#define VTK_COLOR_MODE_DIRECT_SCALARS
int vtkIdType
Definition vtkType.h:354
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:309
#define VTK_SIZEHINT(...)
#define VTK_MARSHALAUTO