source: vtkvis/trunk/Volume.cpp @ 5214

Last change on this file since 5214 was 4650, checked in by ldelgass, 10 years ago

Another fix for OpenGL2 backend

  • Property svn:eol-style set to native
File size: 12.5 KB
Line 
1/* -*- mode: c++; c-basic-offset: 4; indent-tabs-mode: nil -*- */
2/*
3 * Copyright (C) 2004-2012  HUBzero Foundation, LLC
4 *
5 * Author: Leif Delgass <ldelgass@purdue.edu>
6 */
7
8#include <cassert>
9
10#include <vtkDataSet.h>
11#include <vtkPointData.h>
12#include <vtkImageData.h>
13#include <vtkProbeFilter.h>
14#include <vtkVolumeProperty.h>
15#include <vtkGPUVolumeRayCastMapper.h>
16#include <vtkVolumeTextureMapper3D.h>
17#include <vtkRectilinearGrid.h>
18#include <vtkStructuredGrid.h>
19#include <vtkUnstructuredGrid.h>
20#include <vtkPolyData.h>
21#include <vtkCellType.h>
22#include <vtkUnstructuredGridVolumeMapper.h>
23#include <vtkUnstructuredGridVolumeRayCastMapper.h>
24#include <vtkProjectedTetrahedraMapper.h>
25#include <vtkDataSetTriangleFilter.h>
26#include <vtkGaussianSplatter.h>
27
28#include "Volume.h"
29#include "Trace.h"
30
31using namespace VtkVis;
32
33Volume::Volume() :
34    GraphicsObject(),
35    _useUgridMapper(false),
36    _colorMap(NULL)
37{
38}
39
40Volume::~Volume()
41{
42#ifdef WANT_TRACE
43    if (_dataSet != NULL)
44        TRACE("Deleting Volume for %s", _dataSet->getName().c_str());
45    else
46        TRACE("Deleting Volume with NULL DataSet");
47#endif
48}
49
50/**
51 * \brief Create and initialize a VTK Prop to render the Volume
52 */
53void Volume::initProp()
54{
55    if (_prop == NULL) {
56        _prop = vtkSmartPointer<vtkVolume>::New();
57        getVolume()->GetProperty()->SetInterpolationTypeToLinear();
58    }
59}
60
61/**
62 * \brief Get the voxel dimensions (if the volume is not a
63 * uniform grid, unit spacing is returned)
64 */
65void Volume::getSpacing(double spacing[3])
66{
67    spacing[0] = spacing[1] = spacing[2] = 1.0;
68    if (_dataSet == NULL)
69        return;
70    vtkDataSet *ds = _dataSet->getVtkDataSet();
71    if (ds != NULL && vtkImageData::SafeDownCast(ds) != NULL) {
72        vtkImageData::SafeDownCast(ds)->GetSpacing(spacing);
73    } else if (ds != NULL && vtkVolumeMapper::SafeDownCast(_volumeMapper) != NULL) {
74        vtkImageData *imgData = vtkVolumeMapper::SafeDownCast(_volumeMapper)->GetInput();
75        if (imgData != NULL) {
76            imgData->GetSpacing(spacing);
77        }
78    }
79}
80
81/**
82 * \brief Get the average voxel dimension (if the volume is not a
83 * uniform grid, unit spacing is returned)
84 */
85double Volume::getAverageSpacing()
86{
87    double spacing[3];
88    getSpacing(spacing);
89    TRACE("Spacing: %g %g %g", spacing[0], spacing[1], spacing[2]);
90    return (spacing[0] + spacing[1] + spacing[2]) * 0.333;
91}
92
93/**
94 * \brief Internal method to set up pipeline after a state change
95 */
96void Volume::update()
97{
98    if (_dataSet == NULL)
99        return;
100
101    vtkDataSet *ds = _dataSet->getVtkDataSet();
102
103    if (_dataSet->is2D()) {
104        USER_ERROR("Volume rendering requires a 3D data set");
105        _dataSet = NULL;
106        return;
107    }
108
109    if (ds->GetPointData() == NULL ||
110        ds->GetPointData()->GetScalars() == NULL) {
111        USER_ERROR("No scalar field was found in the data set");
112        return;
113    }
114
115    if (vtkImageData::SafeDownCast(ds) != NULL) {
116        // Image data required for these mappers
117#ifdef USE_GPU_RAYCAST_MAPPER
118        _volumeMapper = vtkSmartPointer<vtkGPUVolumeRayCastMapper>::New();
119        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->AutoAdjustSampleDistancesOff();
120#else
121        _volumeMapper = vtkSmartPointer<vtkVolumeTextureMapper3D>::New();
122#endif
123#ifdef USE_VTK6
124#ifdef USE_GPU_RAYCAST_MAPPER
125        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->SetInputData(ds);
126#else
127        vtkVolumeTextureMapper3D::SafeDownCast(_volumeMapper)->SetInputData(ds);
128#endif
129#else
130        _volumeMapper->SetInput(ds);
131#endif
132        vtkVolumeMapper::SafeDownCast(_volumeMapper)->SetBlendModeToComposite();
133    } else if (_dataSet->isCloud()) {
134        // DataSet is a 3D point cloud
135        vtkSmartPointer<vtkGaussianSplatter> splatter = vtkSmartPointer<vtkGaussianSplatter>::New();
136#ifdef USE_VTK6
137        splatter->SetInputData(ds);
138#else
139        splatter->SetInput(ds);
140#endif
141        int dims[3];
142        dims[0] = dims[1] = dims[2] = 64;
143        if (vtkStructuredGrid::SafeDownCast(ds) != NULL) {
144            vtkStructuredGrid::SafeDownCast(ds)->GetDimensions(dims);
145        } else if (vtkRectilinearGrid::SafeDownCast(ds) != NULL) {
146            vtkRectilinearGrid::SafeDownCast(ds)->GetDimensions(dims);
147        }
148        TRACE("Generating volume with dims (%d,%d,%d) from %d points",
149              dims[0], dims[1], dims[2], ds->GetNumberOfPoints());
150        splatter->SetSampleDimensions(dims);
151        splatter->Update();
152        TRACE("Done generating volume");
153#ifdef USE_GPU_RAYCAST_MAPPER
154        _volumeMapper = vtkSmartPointer<vtkGPUVolumeRayCastMapper>::New();
155        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->AutoAdjustSampleDistancesOff();
156#else
157        _volumeMapper = vtkSmartPointer<vtkVolumeTextureMapper3D>::New();
158#endif
159        _volumeMapper->SetInputConnection(splatter->GetOutputPort());
160        vtkVolumeMapper::SafeDownCast(_volumeMapper)->SetBlendModeToComposite();
161    } else if (vtkUnstructuredGrid::SafeDownCast(ds) == NULL || !_useUgridMapper) {
162        // (Slow) Resample using ProbeFilter
163        double bounds[6];
164        ds->GetBounds(bounds);
165        double xLen = bounds[1] - bounds[0];
166        double yLen = bounds[3] - bounds[2];
167        double zLen = bounds[5] - bounds[4];
168
169        int dims[3];
170        dims[0] = dims[1] = dims[2] = 64;
171        if (xLen == 0.0) dims[0] = 1;
172        if (yLen == 0.0) dims[1] = 1;
173        if (zLen == 0.0) dims[2] = 1;
174        if (vtkStructuredGrid::SafeDownCast(ds) != NULL) {
175            vtkStructuredGrid::SafeDownCast(ds)->GetDimensions(dims);
176        } else if (vtkRectilinearGrid::SafeDownCast(ds) != NULL) {
177            vtkRectilinearGrid::SafeDownCast(ds)->GetDimensions(dims);
178        }
179        TRACE("Generating volume with dims (%d,%d,%d) from %d points",
180              dims[0], dims[1], dims[2], ds->GetNumberOfPoints());
181
182        double xSpacing = (dims[0] == 1 ? 0.0 : xLen/((double)(dims[0]-1)));
183        double ySpacing = (dims[1] == 1 ? 0.0 : yLen/((double)(dims[1]-1)));
184        double zSpacing = (dims[2] == 1 ? 0.0 : zLen/((double)(dims[2]-1)));
185
186        vtkSmartPointer<vtkImageData> imageData = vtkSmartPointer<vtkImageData>::New();
187        imageData->SetDimensions(dims[0], dims[1], dims[2]);
188        imageData->SetOrigin(bounds[0], bounds[2], bounds[4]);
189        imageData->SetSpacing(xSpacing, ySpacing, zSpacing);
190
191        vtkSmartPointer<vtkProbeFilter> probe = vtkSmartPointer<vtkProbeFilter>::New();
192        probe->SetInputData(imageData);
193        probe->SetSourceData(ds);
194        probe->Update();
195
196        TRACE("Done generating volume");
197
198#ifdef USE_GPU_RAYCAST_MAPPER
199        _volumeMapper = vtkSmartPointer<vtkGPUVolumeRayCastMapper>::New();
200        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->AutoAdjustSampleDistancesOff();
201#else
202        _volumeMapper = vtkSmartPointer<vtkVolumeTextureMapper3D>::New();
203#endif
204        _volumeMapper->SetInputConnection(probe->GetOutputPort());
205        vtkVolumeMapper::SafeDownCast(_volumeMapper)->SetBlendModeToComposite();
206    } else {
207        // Unstructured grid with cells (not a cloud)
208        vtkUnstructuredGrid *ugrid = vtkUnstructuredGrid::SafeDownCast(ds);
209        assert(ugrid != NULL);
210        // DataSet is unstructured grid
211        // Only works if cells are all tetrahedra
212        //_volumeMapper = vtkSmartPointer<vtkProjectedTetrahedraMapper>::New();
213        // Software raycast rendering - requires all tetrahedra
214        _volumeMapper = vtkSmartPointer<vtkUnstructuredGridVolumeRayCastMapper>::New();
215        if (ugrid->GetCellType(0) == VTK_TETRA &&
216            ugrid->IsHomogeneous()) {
217#ifdef USE_VTK6
218            vtkProjectedTetrahedraMapper::SafeDownCast(_volumeMapper)->SetInputData(ds);
219#else
220            _volumeMapper->SetInput(ds);
221#endif
222        } else {
223            // Decompose to tetrahedra
224            vtkSmartPointer<vtkDataSetTriangleFilter> filter =
225                vtkSmartPointer<vtkDataSetTriangleFilter>::New();
226#ifdef USE_VTK6
227            filter->SetInputData(ugrid);
228#else
229            filter->SetInput(ugrid);
230#endif
231            filter->TetrahedraOnlyOn();
232            TRACE("Decomposing grid to tets");
233            filter->Update();
234            TRACE("Decomposing done");
235            _volumeMapper->SetInputConnection(filter->GetOutputPort());
236        }
237
238        vtkUnstructuredGridVolumeMapper::SafeDownCast(_volumeMapper)->
239            SetBlendModeToComposite();
240    }
241
242    TRACE("Using mapper type: %s", _volumeMapper->GetClassName());
243
244    initProp();
245
246    if (_colorMap == NULL) {
247        setColorMap(ColorMap::getVolumeDefault());
248    }
249
250    setSampleDistance(getAverageSpacing());
251
252    getVolume()->SetMapper(_volumeMapper);
253    _volumeMapper->Update();
254}
255
256void Volume::updateRanges(Renderer *renderer)
257{
258    GraphicsObject::updateRanges(renderer);
259
260    if (getVolume() != NULL) {
261        getVolume()->GetProperty()->SetColor(_colorMap->getColorTransferFunction(_dataRange));
262#ifdef USE_GPU_RAYCAST_MAPPER
263        getVolume()->GetProperty()->SetScalarOpacity(_colorMap->getOpacityTransferFunction(_dataRange, _opacity));
264#else
265        getVolume()->GetProperty()->SetScalarOpacity(_colorMap->getOpacityTransferFunction(_dataRange));
266#endif
267    }
268}
269
270void Volume::updateColorMap()
271{
272    setColorMap(_colorMap);
273}
274
275/**
276 * \brief Assign a color map (transfer function) to use in rendering the Volume
277 */
278void Volume::setColorMap(ColorMap *cmap)
279{
280    if (cmap == NULL)
281        return;
282
283    _colorMap = cmap;
284
285    if (getVolume() != NULL) {
286        getVolume()->GetProperty()->SetColor(_colorMap->getColorTransferFunction(_dataRange));
287#ifdef USE_GPU_RAYCAST_MAPPER
288        getVolume()->GetProperty()->SetScalarOpacity(_colorMap->getOpacityTransferFunction(_dataRange, _opacity));
289#else
290        getVolume()->GetProperty()->SetScalarOpacity(_colorMap->getOpacityTransferFunction(_dataRange));
291#endif
292    }
293}
294
295/**
296 * \brief Return the ColorMap assigned to this Volume
297 */
298ColorMap *Volume::getColorMap()
299{
300    return _colorMap;
301}
302
303/**
304 * \brief Set opacity scaling used to render the Volume
305 */
306void Volume::setOpacity(double opacity)
307{
308    _opacity = opacity;
309    // FIXME: There isn't really a good opacity scaling option that works
310    // across the different mappers/algorithms.  This only works with the
311    // 3D texture mapper, not the GPU raycast mapper
312    if (getVolume() != NULL) {
313#ifdef USE_GPU_RAYCAST_MAPPER
314        getVolume()->GetProperty()->SetScalarOpacity(_colorMap->getOpacityTransferFunction(_dataRange, opacity));
315#else
316        if (opacity < 1.0e-6)
317            opacity = 1.0e-6;
318        getVolume()->GetProperty()->SetScalarOpacityUnitDistance(1.0/opacity);
319#endif
320    }
321}
322
323/**
324 * \brief Set a group of world coordinate planes to clip rendering
325 *
326 * Passing NULL for planes will remove all cliping planes
327 */
328void Volume::setClippingPlanes(vtkPlaneCollection *planes)
329{
330    if (_volumeMapper != NULL) {
331        _volumeMapper->SetClippingPlanes(planes);
332    }
333}
334
335/**
336 * \brief Set distance between samples along rays
337 */
338void Volume::setSampleDistance(float d)
339{
340    TRACE("Sample distance: %g", d);
341    if (_volumeMapper != NULL &&
342#ifdef USE_GPU_RAYCAST_MAPPER
343        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper) != NULL) {
344        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->SetSampleDistance(d);
345#else
346        vtkVolumeTextureMapper3D::SafeDownCast(_volumeMapper) != NULL) {
347        vtkVolumeTextureMapper3D::SafeDownCast(_volumeMapper)->SetSampleDistance(d);
348#endif
349    }
350}
351
352/**
353 * \brief Set Volume renderer blending mode
354 */
355void Volume::setBlendMode(BlendMode mode)
356{
357    if (_volumeMapper == NULL)
358        return;
359
360    vtkVolumeMapper *mapper = vtkVolumeMapper::SafeDownCast(_volumeMapper);
361    if (mapper == NULL) {
362        vtkUnstructuredGridVolumeMapper *ugmapper = vtkUnstructuredGridVolumeMapper::SafeDownCast(_volumeMapper);
363        if (ugmapper == NULL) {
364            TRACE("Mapper does not support BlendMode");
365            return;
366        }
367        switch (mode) {
368        case BLEND_COMPOSITE:
369            ugmapper->SetBlendModeToComposite();
370            break;
371        case BLEND_MAX_INTENSITY:
372            ugmapper->SetBlendModeToMaximumIntensity();
373            break;
374        case BLEND_MIN_INTENSITY:
375        case BLEND_ADDITIVE:
376        default:
377            ERROR("Unknown BlendMode");
378        }
379    } else {
380        switch (mode) {
381        case BLEND_COMPOSITE:
382            mapper->SetBlendModeToComposite();
383            break;
384        case BLEND_MAX_INTENSITY:
385            mapper->SetBlendModeToMaximumIntensity();
386            break;
387        case BLEND_MIN_INTENSITY:
388            mapper->SetBlendModeToMinimumIntensity();
389            break;
390        case BLEND_ADDITIVE:
391            mapper->SetBlendModeToAdditive();
392            break;
393        default:
394            ERROR("Unknown BlendMode");
395        }
396    }
397}
Note: See TracBrowser for help on using the repository browser.