VTK  9.3.0
vtkImageMapper3D.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
29 #ifndef vtkImageMapper3D_h
30 #define vtkImageMapper3D_h
31 
32 #include "vtkAbstractMapper3D.h"
33 #include "vtkRenderingCoreModule.h" // For export macro
34 #include "vtkThreads.h" // for VTK_MAX_THREADS
35 
36 VTK_ABI_NAMESPACE_BEGIN
37 class vtkRenderer;
38 class vtkProp3D;
39 class vtkPoints;
40 class vtkMatrix4x4;
41 class vtkLookupTable;
42 class vtkScalarsToColors;
43 class vtkImageSlice;
44 class vtkImageProperty;
45 class vtkImageData;
46 class vtkMultiThreader;
47 class vtkImageToImageMapper3DFriendship;
48 
49 class VTKRENDERINGCORE_EXPORT vtkImageMapper3D : public vtkAbstractMapper3D
50 {
51 public:
53  void PrintSelf(ostream& os, vtkIndent indent) override;
54 
58  virtual void Render(vtkRenderer* renderer, vtkImageSlice* prop) = 0;
59 
65  void ReleaseGraphicsResources(vtkWindow*) override = 0;
66 
68 
71  void SetInputData(vtkImageData* input);
76 
78 
84  vtkSetMacro(Border, vtkTypeBool);
85  vtkBooleanMacro(Border, vtkTypeBool);
86  vtkGetMacro(Border, vtkTypeBool);
88 
90 
96  vtkSetMacro(Background, vtkTypeBool);
97  vtkBooleanMacro(Background, vtkTypeBool);
98  vtkGetMacro(Background, vtkTypeBool);
100 
102 
107  vtkSetMacro(SliceAtFocalPoint, vtkTypeBool);
108  vtkBooleanMacro(SliceAtFocalPoint, vtkTypeBool);
109  vtkGetMacro(SliceAtFocalPoint, vtkTypeBool);
111 
113 
118  vtkSetMacro(SliceFacesCamera, vtkTypeBool);
119  vtkBooleanMacro(SliceFacesCamera, vtkTypeBool);
120  vtkGetMacro(SliceFacesCamera, vtkTypeBool);
122 
124 
131  vtkGetObjectMacro(SlicePlane, vtkPlane);
133 
139  virtual void GetSlicePlaneInDataCoords(vtkMatrix4x4* propMatrix, double plane[4]);
140 
142 
145  vtkSetClampMacro(NumberOfThreads, int, 1, VTK_MAX_THREADS);
146  vtkGetMacro(NumberOfThreads, int);
148 
150 
159  vtkSetMacro(Streaming, vtkTypeBool);
160  vtkGetMacro(Streaming, vtkTypeBool);
161  vtkBooleanMacro(Streaming, vtkTypeBool);
163 
164  // return the bounds in index space
165  virtual void GetIndexBounds(double extent[6]) = 0;
166 
167 protected:
169  ~vtkImageMapper3D() override;
170 
172 
178 
183  vtkInformation* request, vtkInformationVector** inInfo, vtkInformationVector* outInfo) override;
184 
189  static void CheckerboardRGBA(unsigned char* data, int xsize, int ysize, double originx,
190  double originy, double spacingx, double spacingy);
191 
197  unsigned char* MakeTextureData(vtkImageProperty* property, vtkImageData* input, int extent[6],
198  int& xsize, int& ysize, int& bytesPerPixel, bool& reuseTexture, bool& reuseData);
199 
204  void MakeTextureGeometry(const int extent[6], double coords[12], double tcoords[8]);
205 
213  virtual void ComputeTextureSize(
214  const int extent[6], int& xdim, int& ydim, int imageSize[2], int textureSize[2]);
215 
221 
225  vtkImageSlice* GetCurrentProp() { return this->CurrentProp; }
226 
232 
237  void GetBackgroundColor(vtkImageProperty* property, double color[4]);
238 
245 
246  // The slice.
250 
251  // Information about the image, updated by UpdateInformation
252  double DataSpacing[3];
253  double DataOrigin[3];
254  double DataDirection[9];
255  int DataWholeExtent[6];
256 
257  // Set by vtkImageStack when doing multi-pass rendering
261 
262 private:
263  // The prop this mapper is attached to, or zero if none.
264  vtkImageSlice* CurrentProp;
265  vtkRenderer* CurrentRenderer;
266 
267  // The cached data-to-world matrix
268  vtkMatrix4x4* DataToWorldMatrix;
269 
270  vtkImageMapper3D(const vtkImageMapper3D&) = delete;
271  void operator=(const vtkImageMapper3D&) = delete;
272 
273  friend class vtkImageToImageMapper3DFriendship;
274 };
275 
276 VTK_ABI_NAMESPACE_END
277 #endif
abstract class specifies interface to map 3D data
general representation of visualization data
Definition: vtkDataObject.h:64
abstract class to specify dataset behavior
Definition: vtkDataSet.h:62
topologically and geometrically regular array of data
Definition: vtkImageData.h:52
abstract class for mapping images to the screen
virtual void GetSlicePlaneInDataCoords(vtkMatrix4x4 *propMatrix, double plane[4])
Get the plane as a homogeneous 4-vector that gives the plane equation coefficients.
virtual void ComputeTextureSize(const int extent[6], int &xdim, int &ydim, int imageSize[2], int textureSize[2])
Given an extent that describes a slice (it must have unit thickness in one of the three directions),...
vtkDataObject * GetDataObjectInput()
The input data for this mapper.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static void CheckerboardRGBA(unsigned char *data, int xsize, int ysize, double originx, double originy, double spacingx, double spacingy)
Checkerboard the alpha component of an RGBA image.
vtkScalarsToColors * DefaultLookupTable
vtkMultiThreader * Threader
int FillOutputPortInformation(int port, vtkInformation *info) override
See algorithm for more info.
vtkTypeBool SliceFacesCamera
void GetBackgroundColor(vtkImageProperty *property, double color[4])
Get the background color, by using the first color in the supplied lookup table, or black if there is...
vtkMatrix4x4 * GetDataToWorldMatrix()
Get the data-to-world matrix for this mapper, according to the assembly path for its prop.
virtual void Render(vtkRenderer *renderer, vtkImageSlice *prop)=0
This should only be called by the renderer.
unsigned char * MakeTextureData(vtkImageProperty *property, vtkImageData *input, int extent[6], int &xsize, int &ysize, int &bytesPerPixel, bool &reuseTexture, bool &reuseData)
Perform window/level and color mapping operations to produce unsigned char data that can be used as a...
vtkDataSet * GetDataSetInput()
The input data for this mapper.
vtkTypeBool Background
virtual void GetIndexBounds(double extent[6])=0
void ReleaseGraphicsResources(vtkWindow *) override=0
Release any graphics resources that are being consumed by this mapper.
vtkImageData * GetInput()
The input data for this mapper.
void MakeTextureGeometry(const int extent[6], double coords[12], double tcoords[8])
Compute the coordinates and texture coordinates for the image, given an extent that describes a singl...
~vtkImageMapper3D() override
void SetInputData(vtkImageData *input)
The input data for this mapper.
vtkTypeBool Streaming
vtkTypeBool ProcessRequest(vtkInformation *request, vtkInformationVector **inInfo, vtkInformationVector *outInfo) override
Handle requests from the pipeline executive.
int FillInputPortInformation(int port, vtkInformation *info) override
See algorithm for more info.
vtkRenderer * GetCurrentRenderer()
Get the renderer associated with this mapper, or zero if none.
vtkImageSlice * GetCurrentProp()
Get the vtkImage prop associated with this mapper, or zero if none.
vtkTypeBool SliceAtFocalPoint
image display properties
represents an image in a 3D scene
Definition: vtkImageSlice.h:48
a simple class to control print indentation
Definition: vtkIndent.h:38
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
map scalar values into colors via a lookup table
represent and manipulate 4x4 transformation matrices
Definition: vtkMatrix4x4.h:40
A class for performing multithreaded execution.
perform various plane computations
Definition: vtkPlane.h:35
represent and manipulate 3D points
Definition: vtkPoints.h:38
represents an 3D object for placement in a rendered scene
Definition: vtkProp3D.h:48
abstract specification for renderers
Definition: vtkRenderer.h:71
Superclass for mapping scalar values to colors.
window superclass for vtkRenderWindow
Definition: vtkWindow.h:37
@ Background
Definition: vtkX3D.h:71
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
@ extent
Definition: vtkX3D.h:345
@ color
Definition: vtkX3D.h:221
@ data
Definition: vtkX3D.h:315
int vtkTypeBool
Definition: vtkABI.h:64