VTK  9.2.6
vtkVolumeMapper.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkVolumeMapper.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=========================================================================*/
28#ifndef vtkVolumeMapper_h
29#define vtkVolumeMapper_h
30
32#include "vtkRenderingVolumeModule.h" // For export macro
33
34class vtkImageData;
36class vtkRenderer;
37class vtkVolume;
38
39#define VTK_CROP_SUBVOLUME 0x0002000
40#define VTK_CROP_FENCE 0x2ebfeba
41#define VTK_CROP_INVERTED_FENCE 0x5140145
42#define VTK_CROP_CROSS 0x0417410
43#define VTK_CROP_INVERTED_CROSS 0x7be8bef
44
45class vtkWindow;
46
47class VTKRENDERINGVOLUME_EXPORT vtkVolumeMapper : public vtkAbstractVolumeMapper
48{
49public:
51 void PrintSelf(ostream& os, vtkIndent indent) override;
52
54
57 virtual void SetInputData(vtkImageData*);
58 virtual void SetInputData(vtkDataSet*);
60 virtual vtkDataSet* GetInput();
61 virtual vtkDataSet* GetInput(const int port);
63
65
107 vtkSetMacro(BlendMode, int);
110 {
112 }
114 {
116 }
118 {
120 }
123 void SetBlendModeToSlice() { this->SetBlendMode(vtkVolumeMapper::SLICE_BLEND); }
124 vtkGetMacro(BlendMode, int);
126
128
136 vtkSetVector2Macro(AverageIPScalarRange, double);
137 vtkGetVectorMacro(AverageIPScalarRange, double, 2);
139
141
145 vtkSetClampMacro(Cropping, vtkTypeBool, 0, 1);
146 vtkGetMacro(Cropping, vtkTypeBool);
147 vtkBooleanMacro(Cropping, vtkTypeBool);
149
151
156 vtkSetVector6Macro(CroppingRegionPlanes, double);
157 vtkGetVectorMacro(CroppingRegionPlanes, double, 6);
159
161
165 vtkGetVectorMacro(VoxelCroppingRegionPlanes, double, 6);
167
169
180 vtkSetMacro(ComputeNormalFromOpacity, bool);
181 vtkGetMacro(ComputeNormalFromOpacity, bool);
182 vtkBooleanMacro(ComputeNormalFromOpacity, bool);
184
186
197 vtkSetClampMacro(CroppingRegionFlags, int, 0x0, 0x7ffffff);
198 vtkGetMacro(CroppingRegionFlags, int);
199 void SetCroppingRegionFlagsToSubVolume() { this->SetCroppingRegionFlags(VTK_CROP_SUBVOLUME); }
200 void SetCroppingRegionFlagsToFence() { this->SetCroppingRegionFlags(VTK_CROP_FENCE); }
202 {
203 this->SetCroppingRegionFlags(VTK_CROP_INVERTED_FENCE);
204 }
205 void SetCroppingRegionFlagsToCross() { this->SetCroppingRegionFlags(VTK_CROP_CROSS); }
207 {
208 this->SetCroppingRegionFlags(VTK_CROP_INVERTED_CROSS);
209 }
211
217 void Render(vtkRenderer* ren, vtkVolume* vol) override = 0;
218
226
270 {
277 SLICE_BLEND
278 };
279
280protected:
283
289 double SpacingAdjustedSampleDistance(double inputSpacing[3], int inputExtent[6]);
290
292
296 bool ComputeNormalFromOpacity = false;
297
301 double AverageIPScalarRange[2];
302
304
309 double CroppingRegionPlanes[6];
310 double VoxelCroppingRegionPlanes[6];
314
316
317private:
318 vtkVolumeMapper(const vtkVolumeMapper&) = delete;
319 void operator=(const vtkVolumeMapper&) = delete;
320};
321
322#endif
Abstract class for a volume mapper.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:63
topologically and geometrically regular array of data
Definition: vtkImageData.h:54
a simple class to control print indentation
Definition: vtkIndent.h:40
Store vtkAlgorithm input/output information.
a dataset that is topologically regular with variable spacing in the three coordinate directions
abstract specification for renderers
Definition: vtkRenderer.h:73
Abstract class for a volume mapper.
double SpacingAdjustedSampleDistance(double inputSpacing[3], int inputExtent[6])
Compute a sample distance from the data spacing.
void Render(vtkRenderer *ren, vtkVolume *vol) override=0
WARNING: INTERNAL METHOD - NOT INTENDED FOR GENERAL USE DO NOT USE THIS METHOD OUTSIDE OF THE RENDERI...
void ConvertCroppingRegionPlanesToVoxels()
Cropping variables, and a method for converting the world coordinate cropping region planes to voxel ...
void SetCroppingRegionFlagsToInvertedFence()
Set the flags for the cropping regions.
virtual void SetInputData(vtkDataSet *)
Set/Get the input data.
virtual void SetInputData(vtkImageData *)
Set/Get the input data.
void SetBlendModeToSlice()
Set/Get the blend mode.
virtual vtkDataSet * GetInput(const int port)
Set/Get the input data.
void SetCroppingRegionFlagsToFence()
Set the flags for the cropping regions.
void SetCroppingRegionFlagsToInvertedCross()
Set the flags for the cropping regions.
void SetBlendModeToMinimumIntensity()
Set/Get the blend mode.
void SetBlendModeToAdditive()
Set/Get the blend mode.
int FillInputPortInformation(int, vtkInformation *) override
Fill the input port information objects for this algorithm.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void ReleaseGraphicsResources(vtkWindow *) override
WARNING: INTERNAL METHOD - NOT INTENDED FOR GENERAL USE Release any graphics resources that are being...
void SetBlendModeToMaximumIntensity()
Set/Get the blend mode.
virtual vtkDataSet * GetInput()
Set/Get the input data.
void SetBlendModeToIsoSurface()
Set/Get the blend mode.
void SetCroppingRegionFlagsToCross()
Set the flags for the cropping regions.
BlendModes
Blend modes.
~vtkVolumeMapper() override
void SetBlendModeToComposite()
Set/Get the blend mode.
void SetCroppingRegionFlagsToSubVolume()
Set the flags for the cropping regions.
void SetBlendModeToAverageIntensity()
Set/Get the blend mode.
vtkTypeBool Cropping
Cropping variables, and a method for converting the world coordinate cropping region planes to voxel ...
int CroppingRegionFlags
Cropping variables, and a method for converting the world coordinate cropping region planes to voxel ...
virtual void SetInputData(vtkRectilinearGrid *)
Set/Get the input data.
represents a volume (data & properties) in a rendered scene
Definition: vtkVolume.h:51
window superclass for vtkRenderWindow
Definition: vtkWindow.h:39
int vtkTypeBool
Definition: vtkABI.h:69
#define VTK_CROP_INVERTED_FENCE
#define VTK_CROP_FENCE
#define VTK_CROP_SUBVOLUME
#define VTK_CROP_INVERTED_CROSS
#define VTK_CROP_CROSS