source: trunk/packages/vizservers/vtkvis/RpContour3D.cpp @ 2618

Last change on this file since 2618 was 2612, checked in by ldelgass, 13 years ago

Refactor vtkvis to support setting colormap fields by name/attribute type
rather than always using active scalars/vectors. Also convert common
graphics objects set methods in Renderer to template methods and separate
core and graphics object related methods to separate files.

  • Property svn:eol-style set to native
File size: 6.3 KB
Line 
1/* -*- mode: c++; c-basic-offset: 4; indent-tabs-mode: nil -*- */
2/*
3 * Copyright (C) 2011, Purdue Research Foundation
4 *
5 * Author: Leif Delgass <ldelgass@purdue.edu>
6 */
7
8#include <cassert>
9
10#include <vtkDataSet.h>
11#include <vtkPointData.h>
12#include <vtkCellData.h>
13#include <vtkCellDataToPointData.h>
14#include <vtkContourFilter.h>
15#include <vtkPolyDataMapper.h>
16#include <vtkUnstructuredGrid.h>
17#include <vtkProperty.h>
18#include <vtkDelaunay2D.h>
19#include <vtkDelaunay3D.h>
20#include <vtkDataSetSurfaceFilter.h>
21
22#include "RpContour3D.h"
23#include "Trace.h"
24
25using namespace Rappture::VtkVis;
26
27Contour3D::Contour3D(int numContours) :
28    VtkGraphicsObject(),
29    _numContours(numContours),
30    _colorMap(NULL)
31{
32    _color[0] = 0.0f;
33    _color[1] = 0.0f;
34    _color[2] = 1.0f;
35}
36
37Contour3D::Contour3D(const std::vector<double>& contours) :
38    VtkGraphicsObject(),
39    _numContours(contours.size()),
40    _contours(contours),
41    _colorMap(NULL)
42{
43    _color[0] = 0.0f;
44    _color[1] = 0.0f;
45    _color[2] = 1.0f;
46}
47
48Contour3D::~Contour3D()
49{
50#ifdef WANT_TRACE
51    if (_dataSet != NULL)
52        TRACE("Deleting Contour3D for %s", _dataSet->getName().c_str());
53    else
54        TRACE("Deleting Contour3D with NULL DataSet");
55#endif
56}
57
58/**
59 * \brief Called when the color map has been edited
60 */
61void Contour3D::updateColorMap()
62{
63    setColorMap(_colorMap);
64}
65
66/**
67 * \brief Associate a colormap lookup table with the DataSet
68 */
69void Contour3D::setColorMap(ColorMap *cmap)
70{
71    if (cmap == NULL)
72        return;
73
74    _colorMap = cmap;
75 
76    if (_lut == NULL) {
77        _lut = vtkSmartPointer<vtkLookupTable>::New();
78        if (_contourMapper != NULL) {
79            _contourMapper->UseLookupTableScalarRangeOn();
80            _contourMapper->SetLookupTable(_lut);
81        }
82    }
83
84    _lut->DeepCopy(cmap->getLookupTable());
85    _lut->SetRange(_dataRange);
86}
87
88/**
89 * \brief Internal method to re-compute contours after a state change
90 */
91void Contour3D::update()
92{
93    if (_dataSet == NULL) {
94        return;
95    }
96    vtkDataSet *ds = _dataSet->getVtkDataSet();
97
98    initProp();
99
100    if (_dataSet->is2D()) {
101        ERROR("DataSet is 2D");
102        return;
103    }
104
105    // Contour filter to generate isosurfaces
106    if (_contourFilter == NULL) {
107        _contourFilter = vtkSmartPointer<vtkContourFilter>::New();
108    }
109
110    vtkSmartPointer<vtkCellDataToPointData> cellToPtData;
111
112    if (ds->GetPointData() == NULL ||
113        ds->GetPointData()->GetScalars() == NULL) {
114        WARN("No scalar point data in dataset %s", _dataSet->getName().c_str());
115        if (ds->GetCellData() != NULL &&
116            ds->GetCellData()->GetScalars() != NULL) {
117            cellToPtData =
118                vtkSmartPointer<vtkCellDataToPointData>::New();
119            cellToPtData->SetInput(ds);
120            //cellToPtData->PassCellDataOn();
121            cellToPtData->Update();
122            ds = cellToPtData->GetOutput();
123        } else {
124            ERROR("No scalar cell data in dataset %s", _dataSet->getName().c_str());
125        }
126    }
127
128    vtkPolyData *pd = vtkPolyData::SafeDownCast(ds);
129    if (pd) {
130        // DataSet is a vtkPolyData
131        if (pd->GetNumberOfLines() == 0 &&
132            pd->GetNumberOfPolys() == 0 &&
133            pd->GetNumberOfStrips() == 0) {
134            // DataSet is a point cloud
135            // Generate a 3D unstructured grid
136            vtkSmartPointer<vtkDelaunay3D> mesher = vtkSmartPointer<vtkDelaunay3D>::New();
137            mesher->SetInput(pd);
138            _contourFilter->SetInputConnection(mesher->GetOutputPort());
139        } else {
140            // DataSet is a vtkPolyData with lines and/or polygons
141            ERROR("Not a 3D DataSet");
142            return;
143        }
144    } else {
145         // DataSet is NOT a vtkPolyData
146         _contourFilter->SetInput(ds);
147    }
148
149    _contourFilter->ComputeNormalsOn();
150    _contourFilter->ComputeGradientsOff();
151
152    // Speed up multiple contour computation at cost of extra memory use
153    if (_numContours > 1) {
154        _contourFilter->UseScalarTreeOn();
155    } else {
156        _contourFilter->UseScalarTreeOff();
157    }
158
159    _contourFilter->SetNumberOfContours(_numContours);
160
161    if (_contours.empty()) {
162        // Evenly spaced isovalues
163        _contourFilter->GenerateValues(_numContours, _dataRange[0], _dataRange[1]);
164    } else {
165        // User-supplied isovalues
166        for (int i = 0; i < _numContours; i++) {
167            _contourFilter->SetValue(i, _contours[i]);
168        }
169    }
170    if (_contourMapper == NULL) {
171        _contourMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
172        _contourMapper->SetResolveCoincidentTopologyToPolygonOffset();
173        _contourMapper->SetInputConnection(_contourFilter->GetOutputPort());
174        _contourMapper->SetColorModeToMapScalars();
175        getActor()->SetMapper(_contourMapper);
176    }
177
178    if (_lut == NULL) {
179        setColorMap(ColorMap::getDefault());
180    }
181
182    _contourMapper->Update();
183    TRACE("Contour output %d polys, %d strips",
184          _contourFilter->GetOutput()->GetNumberOfPolys(),
185          _contourFilter->GetOutput()->GetNumberOfStrips());
186}
187
188void Contour3D::updateRanges(Renderer *renderer)
189{
190    VtkGraphicsObject::updateRanges(renderer);
191
192    if (_lut != NULL) {
193        _lut->SetRange(_dataRange);
194    }
195
196    if (_contours.empty() && _numContours > 0) {
197        // Need to recompute isovalues
198        update();
199    }
200}
201
202/**
203 * \brief Specify number of evenly spaced isosurfaces to render
204 *
205 * Will override any existing contours
206 */
207void Contour3D::setContours(int numContours)
208{
209    _contours.clear();
210    _numContours = numContours;
211
212    update();
213}
214
215/**
216 * \brief Specify an explicit list of contour isovalues
217 *
218 * Will override any existing contours
219 */
220void Contour3D::setContourList(const std::vector<double>& contours)
221{
222    _contours = contours;
223    _numContours = (int)_contours.size();
224
225    update();
226}
227
228/**
229 * \brief Get the number of contours
230 */
231int Contour3D::getNumContours() const
232{
233    return _numContours;
234}
235
236/**
237 * \brief Get the contour list (may be empty if number of contours
238 * was specified in place of a list)
239 */
240const std::vector<double>& Contour3D::getContourList() const
241{
242    return _contours;
243}
244
245/**
246 * \brief Set a group of world coordinate planes to clip rendering
247 *
248 * Passing NULL for planes will remove all cliping planes
249 */
250void Contour3D::setClippingPlanes(vtkPlaneCollection *planes)
251{
252    if (_contourMapper != NULL) {
253        _contourMapper->SetClippingPlanes(planes);
254    }
255}
Note: See TracBrowser for help on using the repository browser.