VTK  9.1.0
vtkUnstructuredGridVolumeZSweepMapper.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkUnstructuredGridVolumeZSweepMapper.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=========================================================================*/
35#ifndef vtkUnstructuredGridVolumeZSweepMapper_h
36#define vtkUnstructuredGridVolumeZSweepMapper_h
37
38#include "vtkRenderingVolumeModule.h" // For export macro
40
41class vtkRenderer;
42class vtkVolume;
44class vtkCell;
45class vtkGenericCell;
46class vtkIdList;
48class vtkTransform;
49class vtkMatrix4x4;
51class vtkDoubleArray;
53class vtkRenderWindow;
54
55// Internal classes
57{
58class vtkScreenEdge;
59class vtkSpan;
60class vtkPixelListFrame;
61class vtkUseSet;
62class vtkVertices;
63class vtkSimpleScreenEdge;
64class vtkDoubleScreenEdge;
65class vtkVertexEntry;
66class vtkPixelListEntryMemory;
67};
68
69class VTKRENDERINGVOLUME_EXPORT vtkUnstructuredGridVolumeZSweepMapper
71{
72public:
74 void PrintSelf(ostream& os, vtkIndent indent) override;
75
80
82
87 vtkSetClampMacro(ImageSampleDistance, float, 0.1f, 100.0f);
88 vtkGetMacro(ImageSampleDistance, float);
90
92
96 vtkSetClampMacro(MinimumImageSampleDistance, float, 0.1f, 100.0f);
97 vtkGetMacro(MinimumImageSampleDistance, float);
99
101
105 vtkSetClampMacro(MaximumImageSampleDistance, float, 0.1f, 100.0f);
106 vtkGetMacro(MaximumImageSampleDistance, float);
108
110
116 vtkSetClampMacro(AutoAdjustSampleDistances, vtkTypeBool, 0, 1);
117 vtkGetMacro(AutoAdjustSampleDistances, vtkTypeBool);
118 vtkBooleanMacro(AutoAdjustSampleDistances, vtkTypeBool);
120
122
126 vtkSetClampMacro(IntermixIntersectingGeometry, vtkTypeBool, 0, 1);
127 vtkGetMacro(IntermixIntersectingGeometry, vtkTypeBool);
128 vtkBooleanMacro(IntermixIntersectingGeometry, vtkTypeBool);
130
138
145
147
152 vtkGetObjectMacro(RayIntegrator, vtkUnstructuredGridVolumeRayIntegrator);
154
160 void Render(vtkRenderer* ren, vtkVolume* vol) override;
161
162 vtkGetVectorMacro(ImageInUseSize, int, 2);
163 vtkGetVectorMacro(ImageOrigin, int, 2);
164 vtkGetVectorMacro(ImageViewportSize, int, 2);
165
166protected:
169
174
180
187
192
198
203 void CompositeFunction(double zTarget);
204
208 unsigned char ColorComponentRealToByte(float color);
209
213 void RasterizeFace(vtkIdType faceIds[3], int externalSide);
214
221 void RasterizeTriangle(vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkVertexEntry* ve0,
222 vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkVertexEntry* ve1,
223 vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkVertexEntry* ve2, bool externalFace);
224
231 void RasterizeSpan(int y, vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkScreenEdge* left,
232 vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkScreenEdge* right, bool exitFace);
233
240 void RasterizeLine(vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkVertexEntry* v0,
241 vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkVertexEntry* v1, bool exitFace);
242
243 void StoreRenderTime(vtkRenderer* ren, vtkVolume* vol, float t);
244
246
250 double GetZBufferValue(int x, int y);
251
253
259
265
270
272
277
279
280 // This is how big the image would be if it covered the entire viewport
281 int ImageViewportSize[2];
282
283 // This is how big the allocated memory for image is. This may be bigger
284 // or smaller than ImageFullSize - it will be bigger if necessary to
285 // ensure a power of 2, it will be smaller if the volume only covers a
286 // small region of the viewport
287 int ImageMemorySize[2];
288
289 // This is the size of subregion in ImageSize image that we are using for
290 // the current image. Since ImageSize is a power of 2, there is likely
291 // wasted space in it. This number will be used for things such as clearing
292 // the image if necessary.
293 int ImageInUseSize[2];
294
295 // This is the location in ImageFullSize image where our ImageSize image
296 // is located.
297 int ImageOrigin[2];
298
299 // This is the allocated image
300 unsigned char* Image;
301
302 // This is the accumulating double RGBA image
304
310
312
313 float* ZBuffer;
314 int ZBufferSize[2];
315 int ZBufferOrigin[2];
316
319
320 // if use CellScalars, we need to keep track of the
321 // values on each side of the face and figure out
322 // if the face is used by two cells (twosided) or one cell.
323 double FaceScalars[2];
325
326 vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkSpan* Span;
327 vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkPixelListFrame* PixelListFrame;
328
329 // Used by BuildUseSets().
331
332 vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkUseSet* UseSet;
333
335 vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkVertices* Vertices;
336
339
340 // Used by the main loop
342 int XBounds[2];
343 int YBounds[2];
344
345 vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkSimpleScreenEdge* SimpleEdge;
346 vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkDoubleScreenEdge* DoubleEdge;
347
350
352
353 // Used during compositing
357
358 // Benchmark
360
361 vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkPixelListEntryMemory* MemoryManager;
362
363private:
365 void operator=(const vtkUnstructuredGridVolumeZSweepMapper&) = delete;
366};
367
368#endif
abstract class to specify cell behavior
Definition: vtkCell.h:147
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
dynamic, self-adjusting array of double
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:140
a simple class to control print indentation
Definition: vtkIndent.h:113
represent and manipulate 4x4 transformation matrices
Definition: vtkMatrix4x4.h:145
a list of ids arranged in priority order
helper class that draws the image to the screen
create a window for renderers to draw into
abstract specification for renderers
Definition: vtkRenderer.h:173
record modification and/or execution time
Definition: vtkTimeStamp.h:52
describes linear transformations via a 4x4 matrix
Definition: vtkTransform.h:164
Abstract class for an unstructured grid volume mapper.
a superclass for volume ray integration functions
Unstructured grid volume mapper based the ZSweep Algorithm.
void AllocateUseSet(vtkIdType size)
Allocate an array of usesets of size ‘size’ only if the current one is not large enough.
void RasterizeTriangle(vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkVertexEntry *ve0, vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkVertexEntry *ve1, vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkVertexEntry *ve2, bool externalFace)
Perform scan conversion of a triangle defined by its vertices.
vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkSimpleScreenEdge * SimpleEdge
int ReorderTriangle(vtkIdType v[3], vtkIdType w[3])
Reorder vertices ‘v’ in increasing order in ‘w’.
vtkUnstructuredGridVolumeRayIntegrator * RealRayIntegrator
vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkPixelListFrame * PixelListFrame
float RetrieveRenderTime(vtkRenderer *ren, vtkVolume *vol)
virtual void SetRayIntegrator(vtkUnstructuredGridVolumeRayIntegrator *ri)
Set/Get the helper class for integrating rays.
void RasterizeLine(vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkVertexEntry *v0, vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkVertexEntry *v1, bool exitFace)
Scan conversion of a straight line defined by endpoints v0 and v1.
void MainLoop(vtkRenderWindow *renWin)
MainLoop of the Zsweep algorithm.
unsigned char ColorComponentRealToByte(float color)
Convert and clamp a float color component into a unsigned char.
vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkSpan * Span
void CreateAndCleanPixelList()
Create an empty "pixel list" for each pixel of the screen.
void CompositeFunction(double zTarget)
Do delayed compositing from back to front, stopping at zTarget for each pixel inside the bounding box...
void RasterizeFace(vtkIdType faceIds[3], int externalSide)
Perform scan conversion of a triangle face.
double GetZBufferValue(int x, int y)
Return the value of the z-buffer at screen coordinates (x,y).
void RasterizeSpan(int y, vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkScreenEdge *left, vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkScreenEdge *right, bool exitFace)
Perform scan conversion of an horizontal span from left ro right at line y.
vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkPixelListEntryMemory * MemoryManager
vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkUseSet * UseSet
void AllocateVertices(vtkIdType size)
Allocate a vertex array of size ‘size’ only if the current one is not large enough.
void Render(vtkRenderer *ren, vtkVolume *vol) override
WARNING: INTERNAL METHOD - NOT INTENDED FOR GENERAL USE DO NOT USE THIS METHOD OUTSIDE OF THE RENDERI...
double GetMinimumBoundsDepth(vtkRenderer *ren, vtkVolume *vol)
void SetMaxPixelListSize(int size)
Change the maximum size allowed for a pixel list.
void ProjectAndSortVertices(vtkRenderer *ren, vtkVolume *vol)
Project and sort the vertices by z-coordinates in view space in the "event list" (an heap).
void StoreRenderTime(vtkRenderer *ren, vtkVolume *vol, float t)
void BuildUseSets()
For each vertex, find the list of incident faces.
vtkUnstructuredGridVolumeRayIntegrator * RayIntegrator
vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkVertices * Vertices
static vtkUnstructuredGridVolumeZSweepMapper * New()
Set MaxPixelListSize to 32.
void SavePixelListFrame()
For debugging purpose, save the pixel list frame as a dataset.
int GetMaxPixelListSize()
Maximum size allowed for a pixel list.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkUnstructuredGridVolumeZSweepMapperNamespace::vtkDoubleScreenEdge * DoubleEdge
represents the common properties for rendering a volume.
represents a volume (data & properties) in a rendered scene
Definition: vtkVolume.h:134
@ color
Definition: vtkX3D.h:227
@ size
Definition: vtkX3D.h:259
int vtkTypeBool
Definition: vtkABI.h:69
int vtkIdType
Definition: vtkType.h:332