source: vtkvis/trunk/Contour3D.cpp @ 5835

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

Require VTK >= 6.0.0

  • Property svn:eol-style set to native
File size: 16.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 <vtkVersion.h>
11#include <vtkDataSet.h>
12#include <vtkPointData.h>
13#include <vtkCellData.h>
14#include <vtkTrivialProducer.h>
15#include <vtkCellDataToPointData.h>
16#include <vtkContourFilter.h>
17#include <vtkPolyDataMapper.h>
18#include <vtkUnstructuredGrid.h>
19#include <vtkProperty.h>
20#include <vtkDelaunay2D.h>
21#include <vtkDelaunay3D.h>
22#include <vtkDataSetSurfaceFilter.h>
23
24#include "Contour3D.h"
25#include "Renderer.h"
26#include "Trace.h"
27
28using namespace VtkVis;
29
30Contour3D::Contour3D(int numContours) :
31    GraphicsObject(),
32    _numContours(numContours),
33    _pipelineInitialized(false),
34    _colorMap(NULL),
35    _colorMode(COLOR_BY_SCALAR),
36    _colorFieldType(DataSet::POINT_DATA),
37    _renderer(NULL)
38{
39    _colorFieldRange[0] = DBL_MAX;
40    _colorFieldRange[1] = -DBL_MAX;
41}
42
43Contour3D::Contour3D(const std::vector<double>& contours) :
44    GraphicsObject(),
45    _numContours(contours.size()),
46    _contours(contours),
47    _pipelineInitialized(false),
48    _colorMap(NULL),
49    _colorMode(COLOR_BY_SCALAR),
50    _colorFieldType(DataSet::POINT_DATA),
51    _renderer(NULL)
52{
53    _colorFieldRange[0] = DBL_MAX;
54    _colorFieldRange[1] = -DBL_MAX;
55}
56
57Contour3D::~Contour3D()
58{
59#ifdef WANT_TRACE
60    if (_dataSet != NULL)
61        TRACE("Deleting Contour3D for %s", _dataSet->getName().c_str());
62    else
63        TRACE("Deleting Contour3D with NULL DataSet");
64#endif
65}
66
67void Contour3D::setDataSet(DataSet *dataSet,
68                           Renderer *renderer)
69{
70    if (_dataSet != dataSet) {
71        _dataSet = dataSet;
72        _renderer = renderer;
73
74        if (renderer->getUseCumulativeRange()) {
75            renderer->getCumulativeDataRange(_dataRange,
76                                             _dataSet->getActiveScalarsName(),
77                                             1);
78            const char *activeVectors = _dataSet->getActiveVectorsName();
79            if (activeVectors != NULL) {
80                renderer->getCumulativeDataRange(_vectorMagnitudeRange,
81                                                 _dataSet->getActiveVectorsName(),
82                                                 3);
83                for (int i = 0; i < 3; i++) {
84                    renderer->getCumulativeDataRange(_vectorComponentRange[i],
85                                                     _dataSet->getActiveVectorsName(),
86                                                     3, i);
87                }
88            }
89        } else {
90            _dataSet->getScalarRange(_dataRange);
91            _dataSet->getVectorRange(_vectorMagnitudeRange);
92            for (int i = 0; i < 3; i++) {
93                _dataSet->getVectorRange(_vectorComponentRange[i], i);
94            }
95        }
96
97        update();
98    }
99}
100
101/**
102 * \brief Internal method to re-compute contours after a state change
103 */
104void Contour3D::update()
105{
106    if (_dataSet == NULL) {
107        return;
108    }
109    vtkDataSet *ds = _dataSet->getVtkDataSet();
110
111    if (_dataSet->is2D()) {
112        USER_ERROR("Cannot create isosurface(s) since the data set is not 3D");
113        return;
114    }
115
116    // Contour filter to generate isosurfaces
117    if (_contourFilter == NULL) {
118        _contourFilter = vtkSmartPointer<vtkContourFilter>::New();
119    }
120
121    if (!_pipelineInitialized) {
122        vtkAlgorithmOutput *dsOutput = NULL;
123        vtkSmartPointer<vtkCellDataToPointData> cellToPtData;
124        vtkSmartPointer<vtkTrivialProducer> prod;
125
126        if (ds->GetPointData() == NULL ||
127            ds->GetPointData()->GetScalars() == NULL) {
128            TRACE("No scalar point data in dataset %s", _dataSet->getName().c_str());
129            if (ds->GetCellData() != NULL &&
130                ds->GetCellData()->GetScalars() != NULL) {
131                cellToPtData =
132                    vtkSmartPointer<vtkCellDataToPointData>::New();
133                cellToPtData->SetInputData(ds);
134                cellToPtData->PassCellDataOn();
135                dsOutput = cellToPtData->GetOutputPort();
136            } else {
137                USER_ERROR("No scalar field was found in the data set.");
138            }
139        } else {
140            prod = vtkSmartPointer<vtkTrivialProducer>::New();
141            prod->SetOutput(ds);
142            dsOutput = prod->GetOutputPort();
143        }
144
145        if (_dataSet->isCloud()) {
146            // DataSet is a point cloud
147            // Generate a 3D unstructured grid
148            vtkSmartPointer<vtkDelaunay3D> mesher = vtkSmartPointer<vtkDelaunay3D>::New();
149            mesher->SetInputConnection(dsOutput);
150            _contourFilter->SetInputConnection(mesher->GetOutputPort());
151        } else {
152            // DataSet is 3D with cells.  If DataSet is a surface
153            // (e.g. polydata or ugrid with 2D cells), we will
154            // generate lines instead of surfaces
155            _contourFilter->SetInputConnection(dsOutput);
156        }
157    }
158
159    _pipelineInitialized = true;
160
161    _contourFilter->ComputeNormalsOff();
162    _contourFilter->ComputeGradientsOff();
163
164    // Speed up multiple contour computation at cost of extra memory use
165    if (_numContours > 1) {
166        _contourFilter->UseScalarTreeOn();
167    } else {
168        _contourFilter->UseScalarTreeOff();
169    }
170
171    _contourFilter->SetNumberOfContours(_numContours);
172
173    if (_contours.empty()) {
174        // Evenly spaced isovalues
175        TRACE("Setting %d contours between %g and %g", _numContours, _dataRange[0], _dataRange[1]);
176        _contourFilter->GenerateValues(_numContours, _dataRange[0], _dataRange[1]);
177    } else {
178        // User-supplied isovalues
179        for (int i = 0; i < _numContours; i++) {
180            _contourFilter->SetValue(i, _contours[i]);
181        }
182    }
183
184    initProp();
185
186    if (_normalsGenerator == NULL) {
187        _normalsGenerator = vtkSmartPointer<vtkPolyDataNormals>::New();
188        _normalsGenerator->SetInputConnection(_contourFilter->GetOutputPort());
189        _normalsGenerator->SetFeatureAngle(90.);
190        _normalsGenerator->AutoOrientNormalsOff();
191        _normalsGenerator->ComputePointNormalsOn();
192    }
193
194    if (_mapper == NULL) {
195        _mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
196        _mapper->SetResolveCoincidentTopologyToPolygonOffset();
197        _mapper->SetInputConnection(_normalsGenerator->GetOutputPort());
198        _mapper->ScalarVisibilityOff();
199        //_mapper->SetColorModeToMapScalars();
200        getActor()->SetMapper(_mapper);
201    }
202
203    if (_lut == NULL) {
204        setColorMap(ColorMap::getDefault());
205        setColorMode(_colorMode);
206    }
207
208    _mapper->Update();
209    TRACE("Contour output %d polys, %d strips",
210          _contourFilter->GetOutput()->GetNumberOfPolys(),
211          _contourFilter->GetOutput()->GetNumberOfStrips());
212}
213
214void Contour3D::updateRanges(Renderer *renderer)
215{
216    if (_dataSet == NULL) {
217        ERROR("called before setDataSet");
218        return;
219    }
220
221    if (renderer->getUseCumulativeRange()) {
222        renderer->getCumulativeDataRange(_dataRange,
223                                         _dataSet->getActiveScalarsName(),
224                                         1);
225        const char *activeVectors = _dataSet->getActiveVectorsName();
226        if (activeVectors != NULL) {
227            renderer->getCumulativeDataRange(_vectorMagnitudeRange,
228                                             _dataSet->getActiveVectorsName(),
229                                             3);
230            for (int i = 0; i < 3; i++) {
231                renderer->getCumulativeDataRange(_vectorComponentRange[i],
232                                                 _dataSet->getActiveVectorsName(),
233                                                 3, i);
234            }
235        }
236    } else {
237        _dataSet->getScalarRange(_dataRange);
238        _dataSet->getVectorRange(_vectorMagnitudeRange);
239        for (int i = 0; i < 3; i++) {
240            _dataSet->getVectorRange(_vectorComponentRange[i], i);
241        }
242    }
243 
244    // Need to update color map ranges
245    double *rangePtr = _colorFieldRange;
246    if (_colorFieldRange[0] > _colorFieldRange[1]) {
247        rangePtr = NULL;
248    }
249    setColorMode(_colorMode, _colorFieldType, _colorFieldName.c_str(), rangePtr);
250
251    if (_contours.empty() && _numContours > 0) {
252        // Contour isovalues need to be recomputed
253        update();
254    }
255}
256
257void Contour3D::setContourField(const char *name)
258{
259    if (_contourFilter != NULL) {
260        _contourFilter->SetInputArrayToProcess(0, 0, 0, vtkDataObject::FIELD_ASSOCIATION_POINTS, name);
261        update();
262    }
263}
264
265void Contour3D::setColorMode(ColorMode mode)
266{
267    _colorMode = mode;
268    if (_dataSet == NULL)
269        return;
270
271    switch (mode) {
272    case COLOR_BY_SCALAR:
273        setColorMode(mode,
274                     _dataSet->getActiveScalarsType(),
275                     _dataSet->getActiveScalarsName());
276        break;
277    case COLOR_BY_VECTOR_MAGNITUDE:
278        setColorMode(mode,
279                     _dataSet->getActiveVectorsType(),
280                     _dataSet->getActiveVectorsName());
281        break;
282    case COLOR_BY_VECTOR_X:
283        setColorMode(mode,
284                     _dataSet->getActiveVectorsType(),
285                     _dataSet->getActiveVectorsName());
286        break;
287    case COLOR_BY_VECTOR_Y:
288        setColorMode(mode,
289                     _dataSet->getActiveVectorsType(),
290                     _dataSet->getActiveVectorsName());
291        break;
292    case COLOR_BY_VECTOR_Z:
293        setColorMode(mode,
294                     _dataSet->getActiveVectorsType(),
295                     _dataSet->getActiveVectorsName());
296        break;
297    case COLOR_CONSTANT:
298    default:
299        setColorMode(mode, DataSet::POINT_DATA, NULL, NULL);
300        break;
301    }
302}
303
304void Contour3D::setColorMode(ColorMode mode,
305                             const char *name, double range[2])
306{
307    if (_dataSet == NULL)
308        return;
309    DataSet::DataAttributeType type = DataSet::POINT_DATA;
310    int numComponents = 1;
311    if (name != NULL && strlen(name) > 0 &&
312        !_dataSet->getFieldInfo(name, &type, &numComponents)) {
313        ERROR("Field not found: %s", name);
314        return;
315    }
316    setColorMode(mode, type, name, range);
317}
318
319void Contour3D::setColorMode(ColorMode mode, DataSet::DataAttributeType type,
320                             const char *name, double range[2])
321{
322    _colorMode = mode;
323    _colorFieldType = type;
324    if (name == NULL)
325        _colorFieldName.clear();
326    else
327        _colorFieldName = name;
328    if (range == NULL) {
329        _colorFieldRange[0] = DBL_MAX;
330        _colorFieldRange[1] = -DBL_MAX;
331    } else {
332        memcpy(_colorFieldRange, range, sizeof(double)*2);
333    }
334
335    if (_dataSet == NULL || _mapper == NULL)
336        return;
337
338    switch (type) {
339    case DataSet::POINT_DATA:
340        _mapper->SetScalarModeToUsePointFieldData();
341        break;
342    case DataSet::CELL_DATA:
343        _mapper->SetScalarModeToUseCellFieldData();
344        break;
345    default:
346        ERROR("Unsupported DataAttributeType: %d", type);
347        return;
348    }
349
350    if (name != NULL && strlen(name) > 0) {
351        _mapper->SelectColorArray(name);
352    } else {
353        _mapper->SetScalarModeToDefault();
354    }
355
356    if (_lut != NULL) {
357        if (range != NULL) {
358            _lut->SetRange(range);
359        } else if (name != NULL && strlen(name) > 0) {
360            double r[2];
361            int comp = -1;
362            if (mode == COLOR_BY_VECTOR_X)
363                comp = 0;
364            else if (mode == COLOR_BY_VECTOR_Y)
365                comp = 1;
366            else if (mode == COLOR_BY_VECTOR_Z)
367                comp = 2;
368
369            if (_renderer->getUseCumulativeRange()) {
370                int numComponents;
371                if  (!_dataSet->getFieldInfo(name, type, &numComponents)) {
372                    ERROR("Field not found: %s, type: %d", name, type);
373                    return;
374                } else if (numComponents < comp+1) {
375                    ERROR("Request for component %d in field with %d components",
376                          comp, numComponents);
377                    return;
378                }
379                _renderer->getCumulativeDataRange(r, name, type, numComponents, comp);
380            } else {
381                _dataSet->getDataRange(r, name, type, comp);
382            }
383            _lut->SetRange(r);
384        }
385    }
386
387    switch (mode) {
388    case COLOR_BY_SCALAR:
389        _mapper->ScalarVisibilityOn();
390        break;
391    case COLOR_BY_VECTOR_MAGNITUDE:
392        _mapper->ScalarVisibilityOn();
393        if (_lut != NULL) {
394            _lut->SetVectorModeToMagnitude();
395        }
396        break;
397    case COLOR_BY_VECTOR_X:
398        _mapper->ScalarVisibilityOn();
399        if (_lut != NULL) {
400            _lut->SetVectorModeToComponent();
401            _lut->SetVectorComponent(0);
402        }
403        break;
404    case COLOR_BY_VECTOR_Y:
405        _mapper->ScalarVisibilityOn();
406        if (_lut != NULL) {
407            _lut->SetVectorModeToComponent();
408            _lut->SetVectorComponent(1);
409        }
410        break;
411    case COLOR_BY_VECTOR_Z:
412        _mapper->ScalarVisibilityOn();
413        if (_lut != NULL) {
414            _lut->SetVectorModeToComponent();
415            _lut->SetVectorComponent(2);
416        }
417        break;
418    case COLOR_CONSTANT:
419    default:
420        _mapper->ScalarVisibilityOff();
421        break;
422    }
423}
424
425/**
426 * \brief Called when the color map has been edited
427 */
428void Contour3D::updateColorMap()
429{
430    setColorMap(_colorMap);
431}
432
433/**
434 * \brief Associate a colormap lookup table with the DataSet
435 */
436void Contour3D::setColorMap(ColorMap *cmap)
437{
438    if (cmap == NULL)
439        return;
440
441    _colorMap = cmap;
442 
443    if (_lut == NULL) {
444        _lut = vtkSmartPointer<vtkLookupTable>::New();
445        if (_mapper != NULL) {
446            _mapper->UseLookupTableScalarRangeOn();
447            _mapper->SetLookupTable(_lut);
448        }
449        _lut->DeepCopy(cmap->getLookupTable());
450        switch (_colorMode) {
451        case COLOR_CONSTANT:
452        case COLOR_BY_SCALAR:
453            _lut->SetRange(_dataRange);
454            break;
455        case COLOR_BY_VECTOR_MAGNITUDE:
456            _lut->SetRange(_vectorMagnitudeRange);
457            break;
458        case COLOR_BY_VECTOR_X:
459            _lut->SetRange(_vectorComponentRange[0]);
460            break;
461        case COLOR_BY_VECTOR_Y:
462            _lut->SetRange(_vectorComponentRange[1]);
463            break;
464        case COLOR_BY_VECTOR_Z:
465            _lut->SetRange(_vectorComponentRange[2]);
466            break;
467        default:
468            break;
469        }
470    } else {
471        double range[2];
472        _lut->GetTableRange(range);
473        _lut->DeepCopy(cmap->getLookupTable());
474        _lut->SetRange(range);
475        _lut->Modified();
476    }
477
478    switch (_colorMode) {
479    case COLOR_BY_VECTOR_MAGNITUDE:
480        _lut->SetVectorModeToMagnitude();
481        break;
482    case COLOR_BY_VECTOR_X:
483        _lut->SetVectorModeToComponent();
484        _lut->SetVectorComponent(0);
485        break;
486    case COLOR_BY_VECTOR_Y:
487        _lut->SetVectorModeToComponent();
488        _lut->SetVectorComponent(1);
489        break;
490    case COLOR_BY_VECTOR_Z:
491        _lut->SetVectorModeToComponent();
492        _lut->SetVectorComponent(2);
493        break;
494    default:
495         break;
496    }
497}
498
499/**
500 * \brief Specify number of evenly spaced isosurfaces to render
501 *
502 * Will override any existing contours
503 */
504void Contour3D::setNumContours(int numContours)
505{
506    _contours.clear();
507    _numContours = numContours;
508
509    update();
510}
511
512/**
513 * \brief Specify an explicit list of contour isovalues
514 *
515 * Will override any existing contours
516 */
517void Contour3D::setContourList(const std::vector<double>& contours)
518{
519    _contours = contours;
520    _numContours = (int)_contours.size();
521
522    update();
523}
524
525/**
526 * \brief Get the number of contours
527 */
528int Contour3D::getNumContours() const
529{
530    return _numContours;
531}
532
533/**
534 * \brief Get the contour list (may be empty if number of contours
535 * was specified in place of a list)
536 */
537const std::vector<double>& Contour3D::getContourList() const
538{
539    return _contours;
540}
541
542/**
543 * \brief Set a group of world coordinate planes to clip rendering
544 *
545 * Passing NULL for planes will remove all cliping planes
546 */
547void Contour3D::setClippingPlanes(vtkPlaneCollection *planes)
548{
549    if (_mapper != NULL) {
550        _mapper->SetClippingPlanes(planes);
551    }
552}
Note: See TracBrowser for help on using the repository browser.