source: trunk/packages/vizservers/vtkvis/RpVolume.cpp @ 3596

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

Demote some error about unknown datasets to trace messages

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