source: branches/blt4/packages/vizservers/vtkvis/RpContour2D.cpp @ 2302

Last change on this file since 2302 was 2302, checked in by gah, 13 years ago

update from trunk

File size: 8.6 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 "RpContour2D.h"
23#include "Trace.h"
24
25using namespace Rappture::VtkVis;
26
27Contour2D::Contour2D() :
28    _dataSet(NULL),
29    _numContours(0),
30    _edgeWidth(1.0f),
31    _opacity(1.0)
32{
33    _dataRange[0] = 0;
34    _dataRange[1] = 1;
35    _edgeColor[0] = 0;
36    _edgeColor[1] = 0;
37    _edgeColor[2] = 0;
38}
39
40Contour2D::~Contour2D()
41{
42#ifdef WANT_TRACE
43    if (_dataSet != NULL)
44        TRACE("Deleting Contour2D for %s", _dataSet->getName().c_str());
45    else
46        TRACE("Deleting Contour2D with NULL DataSet");
47#endif
48}
49
50/**
51 * \brief Specify input DataSet used to extract contours
52 */
53void Contour2D::setDataSet(DataSet *dataSet)
54{
55    if (_dataSet != dataSet) {
56        _dataSet = dataSet;
57
58        if (_dataSet != NULL) {
59            double dataRange[2];
60            _dataSet->getDataRange(dataRange);
61            _dataRange[0] = dataRange[0];
62            _dataRange[1] = dataRange[1];
63        }
64
65        update();
66    }
67}
68
69/**
70 * \brief Returns the DataSet this Contour2D renders
71 */
72DataSet *Contour2D::getDataSet()
73{
74    return _dataSet;
75}
76
77/**
78 * \brief Get the VTK Prop for the contour lines
79 */
80vtkProp *Contour2D::getProp()
81{
82    return _contourActor;
83}
84
85/**
86 * \brief Create and initialize a VTK Prop to render isolines
87 */
88void Contour2D::initProp()
89{
90    if (_contourActor == NULL) {
91        _contourActor = vtkSmartPointer<vtkActor>::New();
92        _contourActor->GetProperty()->EdgeVisibilityOn();
93        _contourActor->GetProperty()->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
94        _contourActor->GetProperty()->SetLineWidth(_edgeWidth);
95        _contourActor->GetProperty()->SetOpacity(_opacity);
96        _contourActor->GetProperty()->LightingOff();
97    }
98}
99
100/**
101 * \brief Internal method to re-compute contours after a state change
102 */
103void Contour2D::update()
104{
105    if (_dataSet == NULL) {
106        return;
107    }
108    vtkDataSet *ds = _dataSet->getVtkDataSet();
109
110    initProp();
111
112    // Contour filter to generate isolines
113    if (_contourFilter == NULL) {
114        _contourFilter = vtkSmartPointer<vtkContourFilter>::New();
115    }
116
117    vtkSmartPointer<vtkCellDataToPointData> cellToPtData;
118
119    if (ds->GetPointData() == NULL ||
120        ds->GetPointData()->GetScalars() == NULL) {
121        WARN("No scalar point data in dataset %s", _dataSet->getName().c_str());
122        if (ds->GetCellData() != NULL &&
123            ds->GetCellData()->GetScalars() != NULL) {
124            cellToPtData =
125                vtkSmartPointer<vtkCellDataToPointData>::New();
126            cellToPtData->SetInput(ds);
127            //cellToPtData->PassCellDataOn();
128            cellToPtData->Update();
129            ds = cellToPtData->GetOutput();
130        } else {
131            ERROR("No scalar cell data in dataset %s", _dataSet->getName().c_str());
132        }
133    }
134
135    vtkPolyData *pd = vtkPolyData::SafeDownCast(ds);
136    if (pd) {
137        // DataSet is a vtkPolyData
138        if (pd->GetNumberOfLines() == 0 &&
139            pd->GetNumberOfPolys() == 0 &&
140            pd->GetNumberOfStrips() == 0) {
141            // DataSet is a point cloud
142            if (_dataSet->is2D()) {
143                vtkSmartPointer<vtkDelaunay2D> mesher = vtkSmartPointer<vtkDelaunay2D>::New();
144                mesher->SetInput(pd);
145                _contourFilter->SetInputConnection(mesher->GetOutputPort());
146            } else {
147                vtkSmartPointer<vtkDelaunay3D> mesher = vtkSmartPointer<vtkDelaunay3D>::New();
148                mesher->SetInput(pd);
149                vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
150                gf->SetInputConnection(mesher->GetOutputPort());
151                _contourFilter->SetInputConnection(gf->GetOutputPort());
152            }
153        } else {
154            // DataSet is a vtkPolyData with lines and/or polygons
155            _contourFilter->SetInput(ds);
156        }
157    } else {
158        TRACE("Generating surface for data set");
159        // DataSet is NOT a vtkPolyData
160        vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
161        gf->SetInput(ds);
162        _contourFilter->SetInputConnection(gf->GetOutputPort());
163    }
164
165    _contourFilter->ComputeNormalsOff();
166    _contourFilter->ComputeGradientsOff();
167
168    // Speed up multiple contour computation at cost of extra memory use
169    if (_numContours > 1) {
170        _contourFilter->UseScalarTreeOn();
171    } else {
172        _contourFilter->UseScalarTreeOff();
173    }
174
175    _contourFilter->SetNumberOfContours(_numContours);
176
177    if (_contours.empty()) {
178        // Evenly spaced isovalues
179        _contourFilter->GenerateValues(_numContours, _dataRange[0], _dataRange[1]);
180    } else {
181        // User-supplied isovalues
182        for (int i = 0; i < _numContours; i++) {
183            _contourFilter->SetValue(i, _contours[i]);
184        }
185    }
186    if (_contourMapper == NULL) {
187        _contourMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
188        _contourMapper->SetResolveCoincidentTopologyToPolygonOffset();
189        _contourMapper->SetInputConnection(_contourFilter->GetOutputPort());
190        _contourActor->SetMapper(_contourMapper);
191    }
192
193    _contourMapper->Update();
194}
195
196/**
197 * \brief Specify number of evenly spaced contour lines to render
198 *
199 * Will override any existing contours
200 */
201void Contour2D::setContours(int numContours)
202{
203    _contours.clear();
204    _numContours = numContours;
205
206    if (_dataSet != NULL) {
207        double dataRange[2];
208        _dataSet->getDataRange(dataRange);
209        _dataRange[0] = dataRange[0];
210        _dataRange[1] = dataRange[1];
211    }
212
213    update();
214}
215
216/**
217 * \brief Specify number of evenly spaced contour lines to render
218 * between the given range (including range endpoints)
219 *
220 * Will override any existing contours
221 */
222void Contour2D::setContours(int numContours, double range[2])
223{
224    _contours.clear();
225    _numContours = numContours;
226
227    _dataRange[0] = range[0];
228    _dataRange[1] = range[1];
229
230    update();
231}
232
233/**
234 * \brief Specify an explicit list of contour isovalues
235 *
236 * Will override any existing contours
237 */
238void Contour2D::setContourList(const std::vector<double>& contours)
239{
240    _contours = contours;
241    _numContours = (int)_contours.size();
242
243    update();
244}
245
246/**
247 * \brief Get the number of contours
248 */
249int Contour2D::getNumContours() const
250{
251    return _numContours;
252}
253
254/**
255 * \brief Get the contour list (may be empty if number of contours
256 * was specified in place of a list)
257 */
258const std::vector<double>& Contour2D::getContourList() const
259{
260    return _contours;
261}
262
263/**
264 * \brief Turn on/off rendering of this contour set
265 */
266void Contour2D::setVisibility(bool state)
267{
268    if (_contourActor != NULL) {
269        _contourActor->SetVisibility((state ? 1 : 0));
270    }
271}
272
273/**
274 * \brief Get visibility state of the contour set
275 *
276 * \return Is contour set visible?
277 */
278bool Contour2D::getVisibility() const
279{
280    if (_contourActor == NULL) {
281        return false;
282    } else {
283        return (_contourActor->GetVisibility() != 0);
284    }
285}
286
287/**
288 * \brief Set opacity used to render contour lines
289 */
290void Contour2D::setOpacity(double opacity)
291{
292    _opacity = opacity;
293    if (_contourActor != NULL)
294        _contourActor->GetProperty()->SetOpacity(opacity);
295}
296
297/**
298 * \brief Set RGB color of contour lines
299 */
300void Contour2D::setEdgeColor(float color[3])
301{
302    _edgeColor[0] = color[0];
303    _edgeColor[1] = color[1];
304    _edgeColor[2] = color[2];
305    if (_contourActor != NULL)
306        _contourActor->GetProperty()->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
307}
308
309/**
310 * \brief Set pixel width of contour lines (may be a no-op)
311 */
312void Contour2D::setEdgeWidth(float edgeWidth)
313{
314    _edgeWidth = edgeWidth;
315    if (_contourActor != NULL)
316        _contourActor->GetProperty()->SetLineWidth(_edgeWidth);
317}
318
319/**
320 * \brief Set a group of world coordinate planes to clip rendering
321 *
322 * Passing NULL for planes will remove all cliping planes
323 */
324void Contour2D::setClippingPlanes(vtkPlaneCollection *planes)
325{
326    if (_contourMapper != NULL) {
327        _contourMapper->SetClippingPlanes(planes);
328    }
329}
330
331/**
332 * \brief Turn on/off lighting of this object
333 */
334void Contour2D::setLighting(bool state)
335{
336    if (_contourActor != NULL)
337        _contourActor->GetProperty()->SetLighting((state ? 1 : 0));
338}
Note: See TracBrowser for help on using the repository browser.