source: trunk/packages/vizservers/vtkvis/Volume.cpp @ 3866

Last change on this file since 3866 was 3866, checked in by ldelgass, 11 years ago

Handle rectilinear/structured grids in vtk volume rendering using splat
resampling (for now).

  • Property svn:eol-style set to native
File size: 8.8 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 <vtkVolumeProperty.h>
14#include <vtkGPUVolumeRayCastMapper.h>
15#include <vtkVolumeTextureMapper3D.h>
16#include <vtkRectilinearGrid.h>
17#include <vtkStructuredGrid.h>
18#include <vtkUnstructuredGrid.h>
19#include <vtkPolyData.h>
20#include <vtkCellType.h>
21#include <vtkUnstructuredGridVolumeMapper.h>
22#include <vtkUnstructuredGridVolumeRayCastMapper.h>
23#include <vtkProjectedTetrahedraMapper.h>
24#include <vtkDataSetTriangleFilter.h>
25#include <vtkGaussianSplatter.h>
26
27#include "Volume.h"
28#include "Trace.h"
29
30using namespace VtkVis;
31
32Volume::Volume() :
33    GraphicsObject(),
34    _colorMap(NULL)
35{
36}
37
38Volume::~Volume()
39{
40#ifdef WANT_TRACE
41    if (_dataSet != NULL)
42        TRACE("Deleting Volume for %s", _dataSet->getName().c_str());
43    else
44        TRACE("Deleting Volume with NULL DataSet");
45#endif
46}
47
48/**
49 * \brief Create and initialize a VTK Prop to render the Volume
50 */
51void Volume::initProp()
52{
53    if (_prop == NULL) {
54        _prop = vtkSmartPointer<vtkVolume>::New();
55        getVolume()->GetProperty()->SetInterpolationTypeToLinear();
56    }
57}
58
59/**
60 * \brief Get the voxel dimensions (if the volume is not a
61 * uniform grid, unit spacing is returned)
62 */
63void Volume::getSpacing(double spacing[3])
64{
65    spacing[0] = spacing[1] = spacing[2] = 1.0;
66    if (_dataSet == NULL)
67        return;
68    vtkDataSet *ds = _dataSet->getVtkDataSet();
69    if (ds != NULL && vtkImageData::SafeDownCast(ds) != NULL) {
70        vtkImageData::SafeDownCast(ds)->GetSpacing(spacing);
71    } else if (ds != NULL && vtkVolumeMapper::SafeDownCast(_volumeMapper) != NULL) {
72        vtkImageData *imgData = vtkVolumeMapper::SafeDownCast(_volumeMapper)->GetInput();
73        if (imgData != NULL) {
74            imgData->GetSpacing(spacing);
75        }
76    }
77}
78
79/**
80 * \brief Get the average voxel dimension (if the volume is not a
81 * uniform grid, unit spacing is returned)
82 */
83double Volume::getAverageSpacing()
84{
85    double spacing[3];
86    getSpacing(spacing);
87    TRACE("Spacing: %g %g %g", spacing[0], spacing[1], spacing[2]);
88    return (spacing[0] + spacing[1] + spacing[2]) * 0.333;
89}
90
91/**
92 * \brief Internal method to set up pipeline after a state change
93 */
94void Volume::update()
95{
96    if (_dataSet == NULL)
97        return;
98
99    vtkDataSet *ds = _dataSet->getVtkDataSet();
100
101    if (_dataSet->is2D()) {
102        USER_ERROR("Volume rendering requires a 3D data set");
103        _dataSet = NULL;
104        return;
105    }
106
107    if (vtkImageData::SafeDownCast(ds) != NULL) {
108        // Image data required for these mappers
109#ifdef USE_GPU_RAYCAST_MAPPER
110        _volumeMapper = vtkSmartPointer<vtkGPUVolumeRayCastMapper>::New();
111        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->AutoAdjustSampleDistancesOff();
112#else
113        _volumeMapper = vtkSmartPointer<vtkVolumeTextureMapper3D>::New();
114#endif
115#ifdef USE_VTK6
116#ifdef USE_GPU_RAYCAST_MAPPER
117        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->SetInputData(ds);
118#else
119        vtkVolumeTextureMapper3D::SafeDownCast(_volumeMapper)->SetInputData(ds);
120#endif
121#else
122        _volumeMapper->SetInput(ds);
123#endif
124        vtkVolumeMapper::SafeDownCast(_volumeMapper)->SetBlendModeToComposite();
125    } else if (_dataSet->isCloud() ||
126               vtkUnstructuredGrid::SafeDownCast(ds) == NULL) {
127        // DataSet is a 3D point cloud, rectilinear grid or structured grid
128        vtkSmartPointer<vtkGaussianSplatter> splatter = vtkSmartPointer<vtkGaussianSplatter>::New();
129#ifdef USE_VTK6
130        splatter->SetInputData(ds);
131#else
132        splatter->SetInput(ds);
133#endif
134        int dims[3];
135        dims[0] = dims[1] = dims[2] = 64;
136        if (vtkStructuredGrid::SafeDownCast(ds) != NULL) {
137            vtkStructuredGrid::SafeDownCast(ds)->GetDimensions(dims);
138        } else if (vtkRectilinearGrid::SafeDownCast(ds) != NULL) {
139            vtkRectilinearGrid::SafeDownCast(ds)->GetDimensions(dims);
140        }
141        TRACE("Generating volume with dims (%d,%d,%d) from %d points",
142              dims[0], dims[1], dims[2], ds->GetNumberOfPoints());
143        splatter->SetSampleDimensions(dims);
144        splatter->Update();
145        TRACE("Done generating volume");
146#ifdef USE_GPU_RAYCAST_MAPPER
147        _volumeMapper = vtkSmartPointer<vtkGPUVolumeRayCastMapper>::New();
148        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->AutoAdjustSampleDistancesOff();
149#else
150        _volumeMapper = vtkSmartPointer<vtkVolumeTextureMapper3D>::New();
151#endif
152        _volumeMapper->SetInputConnection(splatter->GetOutputPort());
153        vtkVolumeMapper::SafeDownCast(_volumeMapper)->SetBlendModeToComposite();
154    } else {
155        // Unstructured grid with cells (not a cloud)
156        vtkUnstructuredGrid *ugrid = vtkUnstructuredGrid::SafeDownCast(ds);
157        assert(ugrid != NULL);
158        // DataSet is unstructured grid
159        // Only works if cells are all tetrahedra
160        _volumeMapper = vtkSmartPointer<vtkProjectedTetrahedraMapper>::New();
161        // Software raycast rendering - requires all tetrahedra
162        //_volumeMapper = vtkSmartPointer<vtkUnstructuredGridVolumeRayCastMapper>::New();
163        if (ugrid->GetCellType(0) == VTK_TETRA &&
164            ugrid->IsHomogeneous()) {
165#ifdef USE_VTK6
166            vtkProjectedTetrahedraMapper::SafeDownCast(_volumeMapper)->SetInputData(ds);
167#else
168            _volumeMapper->SetInput(ds);
169#endif
170        } else {
171            // Decompose to tetrahedra
172            vtkSmartPointer<vtkDataSetTriangleFilter> filter =
173                vtkSmartPointer<vtkDataSetTriangleFilter>::New();
174#ifdef USE_VTK6
175            filter->SetInputData(ugrid);
176#else
177            filter->SetInput(ugrid);
178#endif
179            filter->TetrahedraOnlyOn();
180            TRACE("Decomposing grid to tets");
181            filter->Update();
182            TRACE("Decomposing done");
183            _volumeMapper->SetInputConnection(filter->GetOutputPort());
184        }
185
186        vtkUnstructuredGridVolumeMapper::SafeDownCast(_volumeMapper)->
187            SetBlendModeToComposite();
188    }
189
190    TRACE("Using mapper type: %s", _volumeMapper->GetClassName());
191
192    initProp();
193
194    if (ds->GetPointData() == NULL ||
195        ds->GetPointData()->GetScalars() == NULL) {
196        WARN("No scalar point data in dataset %s", _dataSet->getName().c_str());
197    }
198
199    if (_colorMap == NULL) {
200        setColorMap(ColorMap::getVolumeDefault());
201    }
202
203    setSampleDistance(getAverageSpacing());
204
205    getVolume()->SetMapper(_volumeMapper);
206    _volumeMapper->Update();
207}
208
209void Volume::updateRanges(Renderer *renderer)
210{
211    GraphicsObject::updateRanges(renderer);
212
213    if (getVolume() != NULL) {
214        getVolume()->GetProperty()->SetColor(_colorMap->getColorTransferFunction(_dataRange));
215        getVolume()->GetProperty()->SetScalarOpacity(_colorMap->getOpacityTransferFunction(_dataRange));
216    }
217}
218
219void Volume::updateColorMap()
220{
221    setColorMap(_colorMap);
222}
223
224/**
225 * \brief Assign a color map (transfer function) to use in rendering the Volume
226 */
227void Volume::setColorMap(ColorMap *cmap)
228{
229    if (cmap == NULL)
230        return;
231
232    _colorMap = cmap;
233
234    if (getVolume() != NULL) {
235        getVolume()->GetProperty()->SetColor(_colorMap->getColorTransferFunction(_dataRange));
236        getVolume()->GetProperty()->SetScalarOpacity(_colorMap->getOpacityTransferFunction(_dataRange));
237    }
238}
239
240/**
241 * \brief Return the ColorMap assigned to this Volume
242 */
243ColorMap *Volume::getColorMap()
244{
245    return _colorMap;
246}
247
248/**
249 * \brief Set opacity scaling used to render the Volume
250 */
251void Volume::setOpacity(double opacity)
252{
253    _opacity = opacity;
254    // FIXME: There isn't really a good opacity scaling option that works
255    // across the different mappers/algorithms.  This only works with the
256    // 3D texture mapper, not the GPU raycast mapper
257    if (getVolume() != NULL) {
258        if (opacity < 1.0e-6)
259            opacity = 1.0e-6;
260        getVolume()->GetProperty()->SetScalarOpacityUnitDistance(1.0/opacity);
261    }
262}
263
264/**
265 * \brief Set a group of world coordinate planes to clip rendering
266 *
267 * Passing NULL for planes will remove all cliping planes
268 */
269void Volume::setClippingPlanes(vtkPlaneCollection *planes)
270{
271    if (_volumeMapper != NULL) {
272        _volumeMapper->SetClippingPlanes(planes);
273    }
274}
275
276/**
277 * \brief Set distance between samples along rays
278 */
279void Volume::setSampleDistance(float d)
280{
281    TRACE("Sample distance: %g", d);
282    if (_volumeMapper != NULL &&
283#ifdef USE_GPU_RAYCAST_MAPPER
284        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper) != NULL) {
285        vtkGPUVolumeRayCastMapper::SafeDownCast(_volumeMapper)->SetSampleDistance(d);
286#else
287        vtkVolumeTextureMapper3D::SafeDownCast(_volumeMapper) != NULL) {
288        vtkVolumeTextureMapper3D::SafeDownCast(_volumeMapper)->SetSampleDistance(d);
289#endif
290    }
291}
Note: See TracBrowser for help on using the repository browser.