VTK  9.4.20251203
vtkSmartVolumeMapper.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
64#ifndef vtkSmartVolumeMapper_h
65#define vtkSmartVolumeMapper_h
66
67#include "vtkImageReslice.h" // for VTK_RESLICE_NEAREST, VTK_RESLICE_CUBIC
68#include "vtkRenderingVolumeOpenGL2Module.h" // For export macro
69#include "vtkVolumeMapper.h"
70#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
71
72VTK_ABI_NAMESPACE_BEGIN
79class vtkRenderWindow;
80class vtkVolume;
83
84class VTKRENDERINGVOLUMEOPENGL2_EXPORT VTK_MARSHALAUTO vtkSmartVolumeMapper : public vtkVolumeMapper
85{
86public:
89 void PrintSelf(ostream& os, vtkIndent indent) override;
90
92
102 vtkSetMacro(FinalColorWindow, float);
104
106
109 vtkGetMacro(FinalColorWindow, float);
111
113
120 vtkSetMacro(FinalColorLevel, float);
122
124
127 vtkGetMacro(FinalColorLevel, float);
129
130 // The possible values for the default and current render mode ivars
131 enum
132 {
133 DefaultRenderMode = 0,
134 RayCastRenderMode = 1,
135 GPURenderMode = 2,
136 OSPRayRenderMode = 3,
137 AnariRenderMode = 4,
138 UndefinedRenderMode = 5,
139 InvalidRenderMode = 6
140 };
141
147
154
161
168
174
180
182
185 vtkGetMacro(RequestedRenderMode, int);
187
193
195
202 vtkSetMacro(MaxMemoryInBytes, vtkIdType);
203 vtkGetMacro(MaxMemoryInBytes, vtkIdType);
205
207
213 vtkSetClampMacro(MaxMemoryFraction, float, 0.1f, 1.0f);
214 vtkGetMacro(MaxMemoryFraction, float);
216
218
222 vtkSetClampMacro(InterpolationMode, int, VTK_RESLICE_NEAREST, VTK_RESLICE_CUBIC);
223 vtkGetMacro(InterpolationMode, int);
228
235 vtkImageData* image, int blend_mode, double viewDirection[3], double viewUp[3]);
236
238
243 vtkSetClampMacro(UseJittering, vtkTypeBool, 0, 1);
244 vtkGetMacro(UseJittering, vtkTypeBool);
245 vtkBooleanMacro(UseJittering, vtkTypeBool);
247
249
255 vtkSetClampMacro(InteractiveUpdateRate, double, 1.0e-10, 1.0e10);
257
259
264 vtkGetMacro(InteractiveUpdateRate, double);
266
268
276 vtkSetClampMacro(InteractiveAdjustSampleDistances, vtkTypeBool, 0, 1);
277 vtkGetMacro(InteractiveAdjustSampleDistances, vtkTypeBool);
278 vtkBooleanMacro(InteractiveAdjustSampleDistances, vtkTypeBool);
280
282
291 vtkSetClampMacro(AutoAdjustSampleDistances, vtkTypeBool, 0, 1);
292 vtkGetMacro(AutoAdjustSampleDistances, vtkTypeBool);
293 vtkBooleanMacro(AutoAdjustSampleDistances, vtkTypeBool);
295
297
304 vtkSetMacro(SampleDistance, float);
305 vtkGetMacro(SampleDistance, float);
307
309
315 vtkSetClampMacro(GlobalIlluminationReach, float, 0.0f, 1.0f);
316 vtkGetMacro(GlobalIlluminationReach, float);
318
320
326 vtkSetClampMacro(VolumetricScatteringBlending, float, 0.0f, 2.0f);
327 vtkGetMacro(VolumetricScatteringBlending, float);
329
334 void Render(vtkRenderer*, vtkVolume*) override;
335
343
345
353 {
354 DISABLED = -1,
355 MAGNITUDE = 0,
356 COMPONENT = 1,
357 };
358
359 void SetVectorMode(int mode);
360 vtkGetMacro(VectorMode, int);
361
362 vtkSetClampMacro(VectorComponent, int, 0, 3);
363 vtkGetMacro(VectorComponent, int);
365
367
370 vtkSetStringMacro(Transfer2DYAxisArray);
371 vtkGetStringMacro(Transfer2DYAxisArray);
373
375
383 {
384 LowResModeDisabled = 0,
385 LowResModeResample = 1,
386 };
387
388 vtkSetMacro(LowResMode, int);
389 vtkGetMacro(LowResMode, int)
391
392protected:
395
402
409
411
417
419
425
430
432
440
442
451
457
459
471
477
483
488 vtkGetObjectMacro(GPUMapper, vtkGPUVolumeRayCastMapper);
489
491
498
505
510
515
519 float GlobalIlluminationReach = 0.0;
520
524 float VolumetricScatteringBlending = 0.0;
525
531
538
547
549
560
562
569
576
577 int LowResMode = LowResModeDisabled;
578
579private:
581
585 void SetupVectorMode(vtkVolume* vol);
591 void ComputeMagnitudeCellData(vtkDataSet* input, vtkDataArray* arr);
592 void ComputeMagnitudePointData(vtkDataSet* input, vtkDataArray* arr);
594
596 void operator=(const vtkSmartVolumeMapper&) = delete;
597
598 vtkOSPRayVolumeInterface* OSPRayMapper;
599 vtkAnariVolumeInterface* AnariMapper;
600};
601
602VTK_ABI_NAMESPACE_END
603#endif
Removes link dependence on optional ANARI module.
abstract superclass for arrays of numeric data
abstract class to specify dataset behavior
Definition vtkDataSet.h:57
Ray casting performed on the GPU.
topologically and geometrically regular array of data
Resamples an image to be larger or smaller.
a simple class to control print indentation
Definition vtkIndent.h:29
Mapper to render volumes defined as vtkMultiBlockDataSet.
Removes link dependence on optional ospray module.
create a window for renderers to draw into
abstract specification for renderers
Definition vtkRenderer.h:63
Adaptive volume mapper.
float MaxMemoryFraction
GPU mapper-specific memory ivars.
void ReleaseGraphicsResources(vtkWindow *) override
WARNING: INTERNAL METHOD - NOT INTENDED FOR GENERAL USE Release any graphics resources that are being...
vtkTypeBool AutoAdjustSampleDistances
Set whether or not the sample distance should be automatically calculated within the internal volume ...
int Initialized
Initialization variables.
int VectorMode
VectorMode is a special rendering mode for 3-component vectors which makes use of GPURayCastMapper's ...
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int RequestedRenderMode
The requested render mode is used to compute the current render mode.
void SetRequestedRenderModeToAnari()
Set the requested render mode to vtkSmartVolumeMapper::AnariRenderMode.
int GetLastUsedRenderMode()
This will return the render mode used during the previous call to Render().
vtkImageResample * GPUResampleFilter
This is the resample filter that may be used if we need to create a low resolution version of the vol...
vtkFixedPointVolumeRayCastMapper * RayCastMapper
The three potential mappers.
vtkImageData * InputDataMagnitude
This filter is used to compute the magnitude of 3-component data.
float SampleDistance
The distance between sample points along the ray.
vtkDataSet * LastInput
Keep a cache of the last input to the mapper so that input data changes can be propagated to the resa...
vtkIdType MaxMemoryInBytes
GPU mapper-specific memory ivars.
void SetVectorMode(int mode)
VectorMode is a special rendering mode for 3-component vectors which makes use of GPURayCastMapper's ...
void Initialize(vtkRenderer *ren, vtkVolume *vol)
The initialize method.
int InitializedBlendMode
We need to keep track of the blend mode we had when we initialized because we need to reinitialize (a...
void SetInterpolationModeToLinear()
Set interpolation mode for downsampling (lowres GPU) (initial value: cubic).
void ConnectFilterInput(vtkImageResample *f)
Connect input of the vtkSmartVolumeMapper to the input of the internal resample filter by doing a sha...
void SetRequestedRenderModeToOSPRay()
Set the requested render mode to vtkSmartVolumeMapper::OSPRayRenderMode.
void CreateCanonicalView(vtkRenderer *ren, vtkVolume *volume, vtkVolume *volume2, vtkImageData *image, int blend_mode, double viewDirection[3], double viewUp[3])
This method can be used to render a representative view of the input data into the supplied image giv...
int CurrentRenderMode
The requested render mode is used to compute the current render mode.
vtkTypeBool UseJittering
Enable / disable stochastic jittering.
~vtkSmartVolumeMapper() override
vtkDataSet * LastFilterInput
Keep a cache of the last input to the mapper so that input data changes can be propagated to the resa...
void SetInterpolationModeToNearestNeighbor()
Set interpolation mode for downsampling (lowres GPU) (initial value: cubic).
void SetRequestedRenderMode(int mode)
Set the requested render mode.
float FinalColorLevel
Window / level ivars.
static vtkSmartVolumeMapper * New()
double InteractiveUpdateRate
If the DesiredUpdateRate of the vtkRenderWindow causing the Render is at or above this value,...
LowResModeType
LowResDisable disables low res mode (default) LowResResample enable low res mode by automatically res...
void ConnectMapperInput(vtkVolumeMapper *m)
Connect input of the vtkSmartVolumeMapper to the input of the internal volume mapper by doing a shall...
int VectorComponent
VectorMode is a special rendering mode for 3-component vectors which makes use of GPURayCastMapper's ...
vtkTypeBool InteractiveAdjustSampleDistances
If the InteractiveAdjustSampleDistances flag is enabled, vtkSmartVolumeMapper interactively sets and ...
int InterpolationMode
Used for downsampling.
vtkTimeStamp SupportStatusCheckTime
Initialization variables.
void SetRequestedRenderModeToRayCast()
Set the requested render mode to vtkSmartVolumeMapper::RayCastRenderMode.
vtkGPUVolumeRayCastMapper * GPUMapper
The three potential mappers.
void SetRequestedRenderModeToGPU()
Set the requested render mode to vtkSmartVolumeMapper::GPURenderMode.
vtkGPUVolumeRayCastMapper * GPULowResMapper
The three potential mappers.
float FinalColorWindow
Window / level ivars.
vtkTimeStamp MagnitudeUploadTime
VectorMode is a special rendering mode for 3-component vectors which makes use of GPURayCastMapper's ...
int GPUSupported
Initialization variables.
void SetRequestedRenderModeToDefault()
Set the requested render mode to vtkSmartVolumeMapper::DefaultRenderMode.
vtkImageMagnitude * ImageMagnitude
This filter is used to compute the magnitude of 3-component data.
int RayCastSupported
Initialization variables.
void ComputeRenderMode(vtkRenderer *ren, vtkVolume *vol)
The method that computes the render mode from the requested render mode based on the support status f...
void Render(vtkRenderer *, vtkVolume *) override
WARNING: INTERNAL METHOD - NOT INTENDED FOR GENERAL USE Initialize rendering for this volume.
char * Transfer2DYAxisArray
Define the array used for the Y axis of transfer 2D.
void SetInterpolationModeToCubic()
Set interpolation mode for downsampling (lowres GPU) (initial value: cubic).
int LowResGPUNecessary
Initialization variables.
VectorModeType
VectorMode is a special rendering mode for 3-component vectors which makes use of GPURayCastMapper's ...
record modification and/or execution time
Abstract class for a volume mapper.
represents the common properties for rendering a volume.
represents a volume (data & properties) in a rendered scene
Definition vtkVolume.h:41
window superclass for vtkRenderWindow
Definition vtkWindow.h:29
int vtkTypeBool
Definition vtkABI.h:64
#define VTK_RESLICE_CUBIC
#define VTK_RESLICE_NEAREST
int vtkIdType
Definition vtkType.h:315
#define VTK_MARSHALAUTO