source: vtkvis/branches/1.7/Volume.cpp @ 4608

Last change on this file since 4608 was 4608, checked in by ldelgass, 9 years ago

Merge r4243 from trunk: resampling non-uniform grids for volumes

  • Property svn:eol-style set to native
File size: 12.0 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        WARN("No scalar point data in dataset %s", _dataSet->getName().c_str());
112    }
113
114    if (vtkImageData::SafeDownCast(ds) != NULL) {
115        // Image data required for these mappers
116#ifdef USE_GPU_RAYCAST_MAPPER
117        _volumeMapper = vtkSmartPointer<vtkGPUVolumeRayCastMapper>::New();
118        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->AutoAdjustSampleDistancesOff();
119#else
120        _volumeMapper = vtkSmartPointer<vtkVolumeTextureMapper3D>::New();
121#endif
122#ifdef USE_VTK6
123#ifdef USE_GPU_RAYCAST_MAPPER
124        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->SetInputData(ds);
125#else
126        vtkVolumeTextureMapper3D::SafeDownCast(_volumeMapper)->SetInputData(ds);
127#endif
128#else
129        _volumeMapper->SetInput(ds);
130#endif
131        vtkVolumeMapper::SafeDownCast(_volumeMapper)->SetBlendModeToComposite();
132    } else if (_dataSet->isCloud()) {
133        // DataSet is a 3D point cloud
134        vtkSmartPointer<vtkGaussianSplatter> splatter = vtkSmartPointer<vtkGaussianSplatter>::New();
135#ifdef USE_VTK6
136        splatter->SetInputData(ds);
137#else
138        splatter->SetInput(ds);
139#endif
140        int dims[3];
141        dims[0] = dims[1] = dims[2] = 64;
142        if (vtkStructuredGrid::SafeDownCast(ds) != NULL) {
143            vtkStructuredGrid::SafeDownCast(ds)->GetDimensions(dims);
144        } else if (vtkRectilinearGrid::SafeDownCast(ds) != NULL) {
145            vtkRectilinearGrid::SafeDownCast(ds)->GetDimensions(dims);
146        }
147        TRACE("Generating volume with dims (%d,%d,%d) from %d points",
148              dims[0], dims[1], dims[2], ds->GetNumberOfPoints());
149        splatter->SetSampleDimensions(dims);
150        splatter->Update();
151        TRACE("Done generating volume");
152#ifdef USE_GPU_RAYCAST_MAPPER
153        _volumeMapper = vtkSmartPointer<vtkGPUVolumeRayCastMapper>::New();
154        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->AutoAdjustSampleDistancesOff();
155#else
156        _volumeMapper = vtkSmartPointer<vtkVolumeTextureMapper3D>::New();
157#endif
158        _volumeMapper->SetInputConnection(splatter->GetOutputPort());
159        vtkVolumeMapper::SafeDownCast(_volumeMapper)->SetBlendModeToComposite();
160    } else if (vtkUnstructuredGrid::SafeDownCast(ds) == NULL || !_useUgridMapper) {
161        // (Slow) Resample using ProbeFilter
162        double bounds[6];
163        ds->GetBounds(bounds);
164        double xLen = bounds[1] - bounds[0];
165        double yLen = bounds[3] - bounds[2];
166        double zLen = bounds[5] - bounds[4];
167
168        int dims[3];
169        dims[0] = dims[1] = dims[2] = 64;
170        if (xLen == 0.0) dims[0] = 1;
171        if (yLen == 0.0) dims[1] = 1;
172        if (zLen == 0.0) dims[2] = 1;
173        if (vtkStructuredGrid::SafeDownCast(ds) != NULL) {
174            vtkStructuredGrid::SafeDownCast(ds)->GetDimensions(dims);
175        } else if (vtkRectilinearGrid::SafeDownCast(ds) != NULL) {
176            vtkRectilinearGrid::SafeDownCast(ds)->GetDimensions(dims);
177        }
178        TRACE("Generating volume with dims (%d,%d,%d) from %d points",
179              dims[0], dims[1], dims[2], ds->GetNumberOfPoints());
180
181        double xSpacing = (dims[0] == 1 ? 0.0 : xLen/((double)(dims[0]-1)));
182        double ySpacing = (dims[1] == 1 ? 0.0 : yLen/((double)(dims[1]-1)));
183        double zSpacing = (dims[2] == 1 ? 0.0 : zLen/((double)(dims[2]-1)));
184
185        vtkSmartPointer<vtkImageData> imageData = vtkSmartPointer<vtkImageData>::New();
186        imageData->SetDimensions(dims[0], dims[1], dims[2]);
187        imageData->SetOrigin(bounds[0], bounds[2], bounds[4]);
188        imageData->SetSpacing(xSpacing, ySpacing, zSpacing);
189
190        vtkSmartPointer<vtkProbeFilter> probe = vtkSmartPointer<vtkProbeFilter>::New();
191        probe->SetInputData(imageData);
192        probe->SetSourceData(ds);
193        probe->Update();
194
195        TRACE("Done generating volume");
196
197#ifdef USE_GPU_RAYCAST_MAPPER
198        _volumeMapper = vtkSmartPointer<vtkGPUVolumeRayCastMapper>::New();
199        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->AutoAdjustSampleDistancesOff();
200#else
201        _volumeMapper = vtkSmartPointer<vtkVolumeTextureMapper3D>::New();
202#endif
203        _volumeMapper->SetInputConnection(probe->GetOutputPort());
204        vtkVolumeMapper::SafeDownCast(_volumeMapper)->SetBlendModeToComposite();
205    } else {
206        // Unstructured grid with cells (not a cloud)
207        vtkUnstructuredGrid *ugrid = vtkUnstructuredGrid::SafeDownCast(ds);
208        assert(ugrid != NULL);
209        // DataSet is unstructured grid
210        // Only works if cells are all tetrahedra
211        _volumeMapper = vtkSmartPointer<vtkProjectedTetrahedraMapper>::New();
212        // Software raycast rendering - requires all tetrahedra
213        //_volumeMapper = vtkSmartPointer<vtkUnstructuredGridVolumeRayCastMapper>::New();
214        if (ugrid->GetCellType(0) == VTK_TETRA &&
215            ugrid->IsHomogeneous()) {
216#ifdef USE_VTK6
217            vtkProjectedTetrahedraMapper::SafeDownCast(_volumeMapper)->SetInputData(ds);
218#else
219            _volumeMapper->SetInput(ds);
220#endif
221        } else {
222            // Decompose to tetrahedra
223            vtkSmartPointer<vtkDataSetTriangleFilter> filter =
224                vtkSmartPointer<vtkDataSetTriangleFilter>::New();
225#ifdef USE_VTK6
226            filter->SetInputData(ugrid);
227#else
228            filter->SetInput(ugrid);
229#endif
230            filter->TetrahedraOnlyOn();
231            TRACE("Decomposing grid to tets");
232            filter->Update();
233            TRACE("Decomposing done");
234            _volumeMapper->SetInputConnection(filter->GetOutputPort());
235        }
236
237        vtkUnstructuredGridVolumeMapper::SafeDownCast(_volumeMapper)->
238            SetBlendModeToComposite();
239    }
240
241    TRACE("Using mapper type: %s", _volumeMapper->GetClassName());
242
243    initProp();
244
245    if (_colorMap == NULL) {
246        setColorMap(ColorMap::getVolumeDefault());
247    }
248
249    setSampleDistance(getAverageSpacing());
250
251    getVolume()->SetMapper(_volumeMapper);
252    _volumeMapper->Update();
253}
254
255void Volume::updateRanges(Renderer *renderer)
256{
257    GraphicsObject::updateRanges(renderer);
258
259    if (getVolume() != NULL) {
260        getVolume()->GetProperty()->SetColor(_colorMap->getColorTransferFunction(_dataRange));
261#ifdef USE_GPU_RAYCAST_MAPPER
262        getVolume()->GetProperty()->SetScalarOpacity(_colorMap->getOpacityTransferFunction(_dataRange, _opacity));
263#else
264        getVolume()->GetProperty()->SetScalarOpacity(_colorMap->getOpacityTransferFunction(_dataRange));
265#endif
266    }
267}
268
269void Volume::updateColorMap()
270{
271    setColorMap(_colorMap);
272}
273
274/**
275 * \brief Assign a color map (transfer function) to use in rendering the Volume
276 */
277void Volume::setColorMap(ColorMap *cmap)
278{
279    if (cmap == NULL)
280        return;
281
282    _colorMap = cmap;
283
284    if (getVolume() != NULL) {
285        getVolume()->GetProperty()->SetColor(_colorMap->getColorTransferFunction(_dataRange));
286#ifdef USE_GPU_RAYCAST_MAPPER
287        getVolume()->GetProperty()->SetScalarOpacity(_colorMap->getOpacityTransferFunction(_dataRange, _opacity));
288#else
289        getVolume()->GetProperty()->SetScalarOpacity(_colorMap->getOpacityTransferFunction(_dataRange));
290#endif
291    }
292}
293
294/**
295 * \brief Return the ColorMap assigned to this Volume
296 */
297ColorMap *Volume::getColorMap()
298{
299    return _colorMap;
300}
301
302/**
303 * \brief Set opacity scaling used to render the Volume
304 */
305void Volume::setOpacity(double opacity)
306{
307    _opacity = opacity;
308    // FIXME: There isn't really a good opacity scaling option that works
309    // across the different mappers/algorithms.  This only works with the
310    // 3D texture mapper, not the GPU raycast mapper
311    if (getVolume() != NULL) {
312#ifdef USE_GPU_RAYCAST_MAPPER
313        getVolume()->GetProperty()->SetScalarOpacity(_colorMap->getOpacityTransferFunction(_dataRange, opacity));
314#else
315        if (opacity < 1.0e-6)
316            opacity = 1.0e-6;
317        getVolume()->GetProperty()->SetScalarOpacityUnitDistance(1.0/opacity);
318#endif
319    }
320}
321
322/**
323 * \brief Set a group of world coordinate planes to clip rendering
324 *
325 * Passing NULL for planes will remove all cliping planes
326 */
327void Volume::setClippingPlanes(vtkPlaneCollection *planes)
328{
329    if (_volumeMapper != NULL) {
330        _volumeMapper->SetClippingPlanes(planes);
331    }
332}
333
334/**
335 * \brief Set distance between samples along rays
336 */
337void Volume::setSampleDistance(float d)
338{
339    TRACE("Sample distance: %g", d);
340    if (_volumeMapper != NULL &&
341#ifdef USE_GPU_RAYCAST_MAPPER
342        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper) != NULL) {
343        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->SetSampleDistance(d);
344#else
345        vtkVolumeTextureMapper3D::SafeDownCast(_volumeMapper) != NULL) {
346        vtkVolumeTextureMapper3D::SafeDownCast(_volumeMapper)->SetSampleDistance(d);
347#endif
348    }
349}
350
351/**
352 * \brief Set Volume renderer blending mode
353 */
354void Volume::setBlendMode(BlendMode mode)
355{
356    if (_volumeMapper != NULL) {
357        vtkVolumeMapper *mapper = vtkVolumeMapper::SafeDownCast(_volumeMapper);
358        if (mapper == NULL) {
359            TRACE("Mapper does not support BlendMode");
360            return;
361        }
362        switch (mode) {
363        case BLEND_COMPOSITE:
364            mapper->SetBlendModeToComposite();
365            break;
366        case BLEND_MAX_INTENSITY:
367            mapper->SetBlendModeToMaximumIntensity();
368            break;
369        case BLEND_MIN_INTENSITY:
370            mapper->SetBlendModeToMinimumIntensity();
371            break;
372        case BLEND_ADDITIVE:
373            mapper->SetBlendModeToAdditive();
374            break;
375        default:
376            ERROR("Unknown BlendMode");
377        }
378    }
379}
Note: See TracBrowser for help on using the repository browser.