source: trunk/packages/vizservers/vtkvis/RpPseudoColor.cpp @ 2864

Last change on this file since 2864 was 2641, checked in by ldelgass, 13 years ago

Add support for setting glyph scaling and color by named field.

  • Property svn:eol-style set to native
File size: 15.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#include <cfloat>
10#include <cstring>
11
12#include <vtkDataSet.h>
13#include <vtkPointData.h>
14#include <vtkCellData.h>
15#include <vtkDataSetMapper.h>
16#include <vtkUnstructuredGrid.h>
17#include <vtkProperty.h>
18#include <vtkImageData.h>
19#include <vtkLookupTable.h>
20#include <vtkTransform.h>
21#include <vtkDelaunay2D.h>
22#include <vtkDelaunay3D.h>
23#include <vtkGaussianSplatter.h>
24#include <vtkExtractVOI.h>
25#include <vtkDataSetSurfaceFilter.h>
26
27#include "RpPseudoColor.h"
28#include "RpVtkRenderer.h"
29#include "Trace.h"
30
31#define MESH_POINT_CLOUDS
32
33using namespace Rappture::VtkVis;
34
35PseudoColor::PseudoColor() :
36    VtkGraphicsObject(),
37    _colorMap(NULL),
38    _colorMode(COLOR_BY_SCALAR),
39    _colorFieldType(DataSet::POINT_DATA),
40    _renderer(NULL)
41{
42    _colorFieldRange[0] = DBL_MAX;
43    _colorFieldRange[1] = -DBL_MAX;
44}
45
46PseudoColor::~PseudoColor()
47{
48#ifdef WANT_TRACE
49    if (_dataSet != NULL)
50        TRACE("Deleting PseudoColor for %s", _dataSet->getName().c_str());
51    else
52        TRACE("Deleting PseudoColor with NULL DataSet");
53#endif
54}
55
56void PseudoColor::setDataSet(DataSet *dataSet,
57                             Renderer *renderer)
58{
59    if (_dataSet != dataSet) {
60        _dataSet = dataSet;
61
62        _renderer = renderer;
63
64        if (renderer->getUseCumulativeRange()) {
65            renderer->getCumulativeDataRange(_dataRange,
66                                             _dataSet->getActiveScalarsName(),
67                                             1);
68            renderer->getCumulativeDataRange(_vectorMagnitudeRange,
69                                             _dataSet->getActiveVectorsName(),
70                                             3);
71            for (int i = 0; i < 3; i++) {
72                renderer->getCumulativeDataRange(_vectorComponentRange[i],
73                                                 _dataSet->getActiveVectorsName(),
74                                                 3, i);
75            }
76        } else {
77            _dataSet->getScalarRange(_dataRange);
78            _dataSet->getVectorRange(_vectorMagnitudeRange);
79            for (int i = 0; i < 3; i++) {
80                _dataSet->getVectorRange(_vectorComponentRange[i], i);
81            }
82        }
83
84        update();
85    }
86}
87
88/**
89 * \brief Internal method to set up pipeline after a state change
90 */
91void PseudoColor::update()
92{
93    if (_dataSet == NULL)
94        return;
95
96    vtkDataSet *ds = _dataSet->getVtkDataSet();
97
98    // Mapper, actor to render color-mapped data set
99    if (_dsMapper == NULL) {
100        _dsMapper = vtkSmartPointer<vtkDataSetMapper>::New();
101        // Map scalars through lookup table regardless of type
102        _dsMapper->SetColorModeToMapScalars();
103        //_dsMapper->InterpolateScalarsBeforeMappingOn();
104    }
105
106    vtkPolyData *pd = vtkPolyData::SafeDownCast(ds);
107    if (pd) {
108        // DataSet is a vtkPolyData
109        if (pd->GetNumberOfLines() == 0 &&
110            pd->GetNumberOfPolys() == 0 &&
111            pd->GetNumberOfStrips() == 0) {
112            // DataSet is a point cloud
113            PrincipalPlane plane;
114            double offset;
115            if (_dataSet->is2D(&plane, &offset)) {
116#ifdef MESH_POINT_CLOUDS
117                vtkSmartPointer<vtkDelaunay2D> mesher = vtkSmartPointer<vtkDelaunay2D>::New();
118                if (plane == PLANE_ZY) {
119                    vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
120                    trans->RotateWXYZ(90, 0, 1, 0);
121                    if (offset != 0.0) {
122                        trans->Translate(-offset, 0, 0);
123                    }
124                    mesher->SetTransform(trans);
125                } else if (plane == PLANE_XZ) {
126                    vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
127                    trans->RotateWXYZ(-90, 1, 0, 0);
128                    if (offset != 0.0) {
129                        trans->Translate(0, -offset, 0);
130                    }
131                    mesher->SetTransform(trans);
132                } else if (offset != 0.0) {
133                    // XY with Z offset
134                    vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
135                    trans->Translate(0, 0, -offset);
136                    mesher->SetTransform(trans);
137                }
138                mesher->SetInput(pd);
139                _dsMapper->SetInputConnection(mesher->GetOutputPort());
140#else
141                vtkSmartPointer<vtkGaussianSplatter> splatter = vtkSmartPointer<vtkGaussianSplatter>::New();
142                vtkSmartPointer<vtkExtractVOI> slicer = vtkSmartPointer<vtkExtractVOI>::New();
143                splatter->SetInput(pd);
144                int dims[3];
145                splatter->GetSampleDimensions(dims);
146                TRACE("Sample dims: %d %d %d", dims[0], dims[1], dims[2]);
147                if (plane == PLANE_ZY) {
148                    dims[0] = 3;
149                    slicer->SetVOI(1, 1, 0, dims[1]-1, 0, dims[1]-1);
150                } else if (plane == PLANE_XZ) {
151                    dims[1] = 3;
152                    slicer->SetVOI(0, dims[0]-1, 1, 1, 0, dims[2]-1);
153                } else {
154                    dims[2] = 3;
155                    slicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, 1, 1);
156                }
157                splatter->SetSampleDimensions(dims);
158                double bounds[6];
159                splatter->Update();
160                splatter->GetModelBounds(bounds);
161                TRACE("Model bounds: %g %g %g %g %g %g",
162                      bounds[0], bounds[1],
163                      bounds[2], bounds[3],
164                      bounds[4], bounds[5]);
165                slicer->SetInputConnection(splatter->GetOutputPort());
166                slicer->SetSampleRate(1, 1, 1);
167                vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
168                gf->UseStripsOn();
169                gf->SetInputConnection(slicer->GetOutputPort());
170                _dsMapper->SetInputConnection(gf->GetOutputPort());
171#endif
172            } else {
173                vtkSmartPointer<vtkDelaunay3D> mesher = vtkSmartPointer<vtkDelaunay3D>::New();
174                mesher->SetInput(pd);
175                vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
176                gf->SetInputConnection(mesher->GetOutputPort());
177                _dsMapper->SetInputConnection(gf->GetOutputPort());
178             }
179        } else {
180            // DataSet is a vtkPolyData with lines and/or polygons
181            _dsMapper->SetInput(ds);
182        }
183    } else {
184        TRACE("Generating surface for data set");
185        // DataSet is NOT a vtkPolyData
186        vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
187        gf->SetInput(ds);
188        _dsMapper->SetInputConnection(gf->GetOutputPort());
189    }
190
191    if (_lut == NULL) {
192        setColorMap(ColorMap::getDefault());
193    }
194
195    setColorMode(_colorMode);
196
197    initProp();
198    getActor()->SetMapper(_dsMapper);
199    _dsMapper->Update();
200}
201
202void PseudoColor::updateRanges(Renderer *renderer)
203{
204    if (_dataSet == NULL) {
205        ERROR("called before setDataSet");
206        return;
207    }
208
209    if (renderer->getUseCumulativeRange()) {
210        renderer->getCumulativeDataRange(_dataRange,
211                                         _dataSet->getActiveScalarsName(),
212                                         1);
213        renderer->getCumulativeDataRange(_vectorMagnitudeRange,
214                                         _dataSet->getActiveVectorsName(),
215                                         3);
216        for (int i = 0; i < 3; i++) {
217            renderer->getCumulativeDataRange(_vectorComponentRange[i],
218                                             _dataSet->getActiveVectorsName(),
219                                             3, i);
220        }
221    } else {
222        _dataSet->getScalarRange(_dataRange);
223        _dataSet->getVectorRange(_vectorMagnitudeRange);
224        for (int i = 0; i < 3; i++) {
225            _dataSet->getVectorRange(_vectorComponentRange[i], i);
226        }
227    }
228
229    // Need to update color map ranges
230    double *rangePtr = _colorFieldRange;
231    if (_colorFieldRange[0] > _colorFieldRange[1]) {
232        rangePtr = NULL;
233    }
234    setColorMode(_colorMode, _colorFieldType, _colorFieldName.c_str(), rangePtr);
235}
236
237void PseudoColor::setColorMode(ColorMode mode)
238{
239    _colorMode = mode;
240    if (_dataSet == NULL)
241        return;
242
243    switch (mode) {
244    case COLOR_BY_SCALAR:
245        setColorMode(mode,
246                     _dataSet->getActiveScalarsType(),
247                     _dataSet->getActiveScalarsName(),
248                     _dataRange);
249        break;
250    case COLOR_BY_VECTOR_MAGNITUDE:
251        setColorMode(mode,
252                     _dataSet->getActiveVectorsType(),
253                     _dataSet->getActiveVectorsName(),
254                     _vectorMagnitudeRange);
255        break;
256    case COLOR_BY_VECTOR_X:
257        setColorMode(mode,
258                     _dataSet->getActiveVectorsType(),
259                     _dataSet->getActiveVectorsName(),
260                     _vectorComponentRange[0]);
261        break;
262    case COLOR_BY_VECTOR_Y:
263        setColorMode(mode,
264                     _dataSet->getActiveVectorsType(),
265                     _dataSet->getActiveVectorsName(),
266                     _vectorComponentRange[1]);
267        break;
268    case COLOR_BY_VECTOR_Z:
269        setColorMode(mode,
270                     _dataSet->getActiveVectorsType(),
271                     _dataSet->getActiveVectorsName(),
272                     _vectorComponentRange[2]);
273        break;
274    case COLOR_CONSTANT:
275    default:
276        setColorMode(mode, DataSet::POINT_DATA, NULL, NULL);
277        break;
278    }
279}
280
281void PseudoColor::setColorMode(ColorMode mode,
282                               const char *name, double range[2])
283{
284    if (_dataSet == NULL)
285        return;
286    DataSet::DataAttributeType type = DataSet::POINT_DATA;
287    int numComponents = 1;
288    if (name != NULL && strlen(name) > 0 &&
289        !_dataSet->getFieldInfo(name, &type, &numComponents)) {
290        ERROR("Field not found: %s", name);
291        return;
292    }
293    setColorMode(mode, type, name, range);
294}
295
296void PseudoColor::setColorMode(ColorMode mode, DataSet::DataAttributeType type,
297                               const char *name, double range[2])
298{
299    _colorMode = mode;
300    _colorFieldType = type;
301    if (name == NULL)
302        _colorFieldName.clear();
303    else
304        _colorFieldName = name;
305    if (range == NULL) {
306        _colorFieldRange[0] = DBL_MAX;
307        _colorFieldRange[1] = -DBL_MAX;
308    } else {
309        memcpy(_colorFieldRange, range, sizeof(double)*2);
310    }
311
312    if (_dataSet == NULL || _dsMapper == NULL)
313        return;
314
315    switch (type) {
316    case DataSet::POINT_DATA:
317        _dsMapper->SetScalarModeToUsePointFieldData();
318        break;
319    case DataSet::CELL_DATA:
320        _dsMapper->SetScalarModeToUseCellFieldData();
321        break;
322    default:
323        ERROR("Unsupported DataAttributeType: %d", type);
324        return;
325    }
326
327    if (name != NULL && strlen(name) > 0) {
328        _dsMapper->SelectColorArray(name);
329    } else {
330        _dsMapper->SetScalarModeToDefault();
331    }
332
333    if (_lut != NULL) {
334        if (range != NULL) {
335            _lut->SetRange(range);
336        } else if (name != NULL && strlen(name) > 0) {
337            double r[2];
338            int comp = -1;
339            if (mode == COLOR_BY_VECTOR_X)
340                comp = 0;
341            else if (mode == COLOR_BY_VECTOR_Y)
342                comp = 1;
343            else if (mode == COLOR_BY_VECTOR_Z)
344                comp = 2;
345
346            if (_renderer->getUseCumulativeRange()) {
347                int numComponents;
348                if  (!_dataSet->getFieldInfo(name, type, &numComponents)) {
349                    ERROR("Field not found: %s, type: %d", name, type);
350                    return;
351                } else if (numComponents < comp+1) {
352                    ERROR("Request for component %d in field with %d components",
353                          comp, numComponents);
354                    return;
355                }
356                _renderer->getCumulativeDataRange(r, name, type, numComponents, comp);
357            } else {
358                _dataSet->getDataRange(r, name, type, comp);
359            }
360            _lut->SetRange(r);
361        }
362    }
363
364    switch (mode) {
365    case COLOR_BY_SCALAR:
366        _dsMapper->ScalarVisibilityOn();
367        break;
368    case COLOR_BY_VECTOR_MAGNITUDE:
369        _dsMapper->ScalarVisibilityOn();
370        if (_lut != NULL) {
371            _lut->SetVectorModeToMagnitude();
372        }
373        break;
374    case COLOR_BY_VECTOR_X:
375        _dsMapper->ScalarVisibilityOn();
376        if (_lut != NULL) {
377            _lut->SetVectorModeToComponent();
378            _lut->SetVectorComponent(0);
379        }
380        break;
381    case COLOR_BY_VECTOR_Y:
382        _dsMapper->ScalarVisibilityOn();
383        if (_lut != NULL) {
384            _lut->SetVectorModeToComponent();
385            _lut->SetVectorComponent(1);
386        }
387        break;
388    case COLOR_BY_VECTOR_Z:
389        _dsMapper->ScalarVisibilityOn();
390        if (_lut != NULL) {
391            _lut->SetVectorModeToComponent();
392            _lut->SetVectorComponent(2);
393        }
394        break;
395    case COLOR_CONSTANT:
396    default:
397        _dsMapper->ScalarVisibilityOff();
398        break;
399    }
400}
401
402/**
403 * \brief Called when the color map has been edited
404 */
405void PseudoColor::updateColorMap()
406{
407    setColorMap(_colorMap);
408}
409
410/**
411 * \brief Associate a colormap lookup table with the DataSet
412 */
413void PseudoColor::setColorMap(ColorMap *cmap)
414{
415    if (cmap == NULL)
416        return;
417
418    _colorMap = cmap;
419 
420    if (_lut == NULL) {
421        _lut = vtkSmartPointer<vtkLookupTable>::New();
422        if (_dsMapper != NULL) {
423            _dsMapper->UseLookupTableScalarRangeOn();
424            _dsMapper->SetLookupTable(_lut);
425        }
426        _lut->DeepCopy(cmap->getLookupTable());
427        switch (_colorMode) {
428        case COLOR_CONSTANT:
429        case COLOR_BY_SCALAR:
430            _lut->SetRange(_dataRange);
431            break;
432        case COLOR_BY_VECTOR_MAGNITUDE:
433            _lut->SetRange(_vectorMagnitudeRange);
434            break;
435        case COLOR_BY_VECTOR_X:
436            _lut->SetRange(_vectorComponentRange[0]);
437            break;
438        case COLOR_BY_VECTOR_Y:
439            _lut->SetRange(_vectorComponentRange[1]);
440            break;
441        case COLOR_BY_VECTOR_Z:
442            _lut->SetRange(_vectorComponentRange[2]);
443            break;
444        default:
445            break;
446        }
447    } else {
448        double range[2];
449        _lut->GetTableRange(range);
450        _lut->DeepCopy(cmap->getLookupTable());
451        _lut->SetRange(range);
452    }
453
454    switch (_colorMode) {
455    case COLOR_BY_VECTOR_MAGNITUDE:
456        _lut->SetVectorModeToMagnitude();
457        break;
458    case COLOR_BY_VECTOR_X:
459        _lut->SetVectorModeToComponent();
460        _lut->SetVectorComponent(0);
461        break;
462    case COLOR_BY_VECTOR_Y:
463        _lut->SetVectorModeToComponent();
464        _lut->SetVectorComponent(1);
465        break;
466    case COLOR_BY_VECTOR_Z:
467        _lut->SetVectorModeToComponent();
468        _lut->SetVectorComponent(2);
469        break;
470    default:
471         break;
472    }
473}
474
475/**
476 * \brief Set a group of world coordinate planes to clip rendering
477 *
478 * Passing NULL for planes will remove all cliping planes
479 */
480void PseudoColor::setClippingPlanes(vtkPlaneCollection *planes)
481{
482    if (_dsMapper != NULL) {
483        _dsMapper->SetClippingPlanes(planes);
484    }
485}
Note: See TracBrowser for help on using the repository browser.