source: vtkvis/trunk/PseudoColor.cpp

Last change on this file was 6458, checked in by ldelgass, 8 years ago

Apply fix from r6408 (setting colormode with empty field name) across all
graphics objects.

  • Property svn:eol-style set to native
File size: 17.8 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#include <cfloat>
10#include <cstring>
11
12#include <vtkDataSet.h>
13#include <vtkPointData.h>
14#include <vtkCellData.h>
15#include <vtkPolyData.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#include <vtkVertexGlyphFilter.h>
27
28#include "PseudoColor.h"
29#include "Renderer.h"
30#include "Trace.h"
31
32using namespace VtkVis;
33
34PseudoColor::PseudoColor() :
35    GraphicsObject(),
36    _colorMap(NULL),
37    _colorMode(COLOR_BY_SCALAR),
38    _colorFieldType(DataSet::POINT_DATA),
39    _renderer(NULL),
40    _cloudStyle(CLOUD_MESH)
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 (_mapper == NULL) {
100        _mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
101        _mapper->SetResolveCoincidentTopologyToPolygonOffset();
102        // Map scalars through lookup table regardless of type
103        _mapper->SetColorModeToMapScalars();
104    }
105
106    _splatter = NULL;
107    if (_dataSet->isCloud()) {
108        // DataSet is a point cloud
109        if (_cloudStyle == CLOUD_POINTS ||
110            _dataSet->numDimensions() < 2 || ds->GetNumberOfPoints() < 3) {
111            vtkSmartPointer<vtkVertexGlyphFilter> vgf = vtkSmartPointer<vtkVertexGlyphFilter>::New();
112            vgf->SetInputData(ds);
113            vgf->ReleaseDataFlagOn();
114            _mapper->SetInputConnection(vgf->GetOutputPort());
115        } else {
116            PrincipalPlane plane;
117            double offset;
118            if (_dataSet->is2D(&plane, &offset)) {
119                if (_cloudStyle == CLOUD_MESH) {
120                    // 2D CLOUD_MESH
121                    vtkSmartPointer<vtkDelaunay2D> mesher = vtkSmartPointer<vtkDelaunay2D>::New();
122                    if (plane == PLANE_ZY) {
123                        vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
124                        trans->RotateWXYZ(90, 0, 1, 0);
125                        if (offset != 0.0) {
126                            trans->Translate(-offset, 0, 0);
127                        }
128                        mesher->SetTransform(trans);
129                    } else if (plane == PLANE_XZ) {
130                        vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
131                        trans->RotateWXYZ(-90, 1, 0, 0);
132                        if (offset != 0.0) {
133                            trans->Translate(0, -offset, 0);
134                        }
135                        mesher->SetTransform(trans);
136                    } else if (offset != 0.0) {
137                        // XY with Z offset
138                        vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
139                        trans->Translate(0, 0, -offset);
140                        mesher->SetTransform(trans);
141                    }
142                    mesher->SetInputData(ds);
143                    _mapper->SetInputConnection(mesher->GetOutputPort());
144                } else {
145                    // 2D CLOUD_SPLAT
146                    if (_splatter == NULL) {
147                        _splatter = vtkSmartPointer<vtkGaussianSplatter>::New();
148                    }
149                    _splatter->SetInputData(ds);
150                    int dims[3];
151                    _splatter->GetSampleDimensions(dims);
152                    TRACE("Sample dims: %d %d %d", dims[0], dims[1], dims[2]);
153                    vtkSmartPointer<vtkExtractVOI> slicer = vtkSmartPointer<vtkExtractVOI>::New();
154                    if (plane == PLANE_ZY) {
155                        dims[0] = 3;
156                        slicer->SetVOI(1, 1, 0, dims[1]-1, 0, dims[1]-1);
157                    } else if (plane == PLANE_XZ) {
158                        dims[1] = 3;
159                        slicer->SetVOI(0, dims[0]-1, 1, 1, 0, dims[2]-1);
160                    } else {
161                        dims[2] = 3;
162                        slicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, 1, 1);
163                    }
164                    _splatter->SetSampleDimensions(dims);
165                    double bounds[6];
166                    _splatter->Update();
167                    _splatter->GetModelBounds(bounds);
168                    TRACE("Model bounds: %g %g %g %g %g %g",
169                          bounds[0], bounds[1],
170                          bounds[2], bounds[3],
171                          bounds[4], bounds[5]);
172                    slicer->SetInputConnection(_splatter->GetOutputPort());
173                    slicer->SetSampleRate(1, 1, 1);
174                    vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
175                    gf->UseStripsOn();
176                    gf->SetInputConnection(slicer->GetOutputPort());
177                    _mapper->SetInputConnection(gf->GetOutputPort());
178                }
179            } else {
180                // 3D cloud
181                vtkSmartPointer<vtkDelaunay3D> mesher = vtkSmartPointer<vtkDelaunay3D>::New();
182                mesher->SetInputData(ds);
183                vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
184                gf->SetInputConnection(mesher->GetOutputPort());
185                _mapper->SetInputConnection(gf->GetOutputPort());
186            }
187        }
188    } else if (vtkPolyData::SafeDownCast(ds) != NULL) {
189        // DataSet is a vtkPolyData with lines and/or polygons
190        _mapper->SetInputData(vtkPolyData::SafeDownCast(ds));
191    } else {
192        TRACE("Generating surface for data set");
193        // DataSet is NOT a vtkPolyData and has cells
194        vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
195        gf->SetInputData(ds);
196        _mapper->SetInputConnection(gf->GetOutputPort());
197    }
198
199    setInterpolateBeforeMapping(true);
200
201    if (_lut == NULL) {
202        setColorMap(ColorMap::getDefault());
203        if (ds->GetPointData()->GetScalars() == NULL &&
204            ds->GetPointData()->GetVectors() != NULL) {
205            TRACE("Setting color mode to vector magnitude");
206            setColorMode(COLOR_BY_VECTOR_MAGNITUDE);
207        } else {
208            TRACE("Setting color mode to scalar");
209            setColorMode(COLOR_BY_SCALAR);
210        }
211    } else {
212        double *rangePtr = _colorFieldRange;
213        if (_colorFieldRange[0] > _colorFieldRange[1]) {
214            rangePtr = NULL;
215        }
216        setColorMode(_colorMode, _colorFieldType, _colorFieldName.c_str(), rangePtr);
217    }
218
219    initProp();
220    getActor()->SetMapper(_mapper);
221    _mapper->Update();
222#ifdef WANT_TRACE
223    double *b = getBounds();
224    TRACE("bounds: %g %g %g %g %g %g", b[0], b[1], b[2], b[3], b[4], b[5]);
225#endif
226}
227
228void PseudoColor::setCloudStyle(CloudStyle style)
229{
230    if (style != _cloudStyle) {
231        _cloudStyle = style;
232        if (_dataSet != NULL) {
233            update();
234        }
235    }
236}
237
238void PseudoColor::setInterpolateBeforeMapping(bool state)
239{
240    if (_mapper != NULL) {
241        _mapper->SetInterpolateScalarsBeforeMapping((state ? 1 : 0));
242    }
243}
244
245void PseudoColor::updateRanges(Renderer *renderer)
246{
247    if (_dataSet == NULL) {
248        ERROR("called before setDataSet");
249        return;
250    }
251
252    if (renderer->getUseCumulativeRange()) {
253        renderer->getCumulativeDataRange(_dataRange,
254                                         _dataSet->getActiveScalarsName(),
255                                         1);
256        renderer->getCumulativeDataRange(_vectorMagnitudeRange,
257                                         _dataSet->getActiveVectorsName(),
258                                         3);
259        for (int i = 0; i < 3; i++) {
260            renderer->getCumulativeDataRange(_vectorComponentRange[i],
261                                             _dataSet->getActiveVectorsName(),
262                                             3, i);
263        }
264    } else {
265        _dataSet->getScalarRange(_dataRange);
266        _dataSet->getVectorRange(_vectorMagnitudeRange);
267        for (int i = 0; i < 3; i++) {
268            _dataSet->getVectorRange(_vectorComponentRange[i], i);
269        }
270    }
271
272    // Need to update color map ranges
273    double *rangePtr = _colorFieldRange;
274    if (_colorFieldRange[0] > _colorFieldRange[1]) {
275        rangePtr = NULL;
276    }
277    setColorMode(_colorMode, _colorFieldType, _colorFieldName.c_str(), rangePtr);
278}
279
280void PseudoColor::setColorMode(ColorMode mode)
281{
282    _colorMode = mode;
283    if (_dataSet == NULL)
284        return;
285
286    switch (mode) {
287    case COLOR_BY_SCALAR:
288        setColorMode(mode,
289                     _dataSet->getActiveScalarsType(),
290                     _dataSet->getActiveScalarsName());
291        break;
292    case COLOR_BY_VECTOR_MAGNITUDE:
293        setColorMode(mode,
294                     _dataSet->getActiveVectorsType(),
295                     _dataSet->getActiveVectorsName());
296        break;
297    case COLOR_BY_VECTOR_X:
298        setColorMode(mode,
299                     _dataSet->getActiveVectorsType(),
300                     _dataSet->getActiveVectorsName());
301        break;
302    case COLOR_BY_VECTOR_Y:
303        setColorMode(mode,
304                     _dataSet->getActiveVectorsType(),
305                     _dataSet->getActiveVectorsName());
306        break;
307    case COLOR_BY_VECTOR_Z:
308        setColorMode(mode,
309                     _dataSet->getActiveVectorsType(),
310                     _dataSet->getActiveVectorsName());
311        break;
312    case COLOR_CONSTANT:
313    default:
314        setColorMode(mode, DataSet::POINT_DATA, NULL, NULL);
315        break;
316    }
317}
318
319void PseudoColor::setColorMode(ColorMode mode,
320                               const char *name, double range[2])
321{
322    if (_dataSet == NULL)
323        return;
324    DataSet::DataAttributeType type = DataSet::POINT_DATA;
325    int numComponents = 1;
326    if (name != NULL && strlen(name) > 0 &&
327        !_dataSet->getFieldInfo(name, &type, &numComponents)) {
328        ERROR("Field not found: %s", name);
329        return;
330    }
331    setColorMode(mode, type, name, range);
332}
333
334void PseudoColor::setColorMode(ColorMode mode, DataSet::DataAttributeType type,
335                               const char *name, double range[2])
336{
337    _colorMode = mode;
338    _colorFieldType = type;
339    if (name == NULL)
340        _colorFieldName.clear();
341    else
342        _colorFieldName = name;
343    if (range == NULL) {
344        _colorFieldRange[0] = DBL_MAX;
345        _colorFieldRange[1] = -DBL_MAX;
346    } else {
347        memcpy(_colorFieldRange, range, sizeof(double)*2);
348    }
349
350    if (_dataSet == NULL || _mapper == NULL)
351        return;
352
353    switch (type) {
354    case DataSet::POINT_DATA:
355        _mapper->SetScalarModeToUsePointFieldData();
356        break;
357    case DataSet::CELL_DATA:
358        _mapper->SetScalarModeToUseCellFieldData();
359        break;
360    case DataSet::FIELD_DATA:
361        _mapper->SetScalarModeToUseFieldData();
362        break;
363    default:
364        ERROR("Unsupported DataAttributeType: %d", type);
365        return;
366    }
367
368    if (_splatter != NULL) {
369        if (name != NULL && strlen(name) > 0) {
370            _splatter->SetInputArrayToProcess(0, 0, 0, type, name);
371        }
372        _mapper->SelectColorArray("SplatterValues");
373    } else if (name != NULL && strlen(name) > 0) {
374        _mapper->SelectColorArray(name);
375    } else {
376        _mapper->SetScalarModeToDefault();
377    }
378
379    if (_lut != NULL) {
380        if (range != NULL) {
381            _lut->SetRange(range);
382        } else if (name != NULL && strlen(name) > 0) {
383            double r[2];
384            int comp = -1;
385            if (mode == COLOR_BY_VECTOR_X)
386                comp = 0;
387            else if (mode == COLOR_BY_VECTOR_Y)
388                comp = 1;
389            else if (mode == COLOR_BY_VECTOR_Z)
390                comp = 2;
391
392            if (_renderer->getUseCumulativeRange()) {
393                int numComponents;
394                if  (!_dataSet->getFieldInfo(name, type, &numComponents)) {
395                    ERROR("Field not found: %s, type: %d", name, type);
396                    return;
397                } else if (numComponents < comp+1) {
398                    ERROR("Request for component %d in field with %d components",
399                          comp, numComponents);
400                    return;
401                }
402                _renderer->getCumulativeDataRange(r, name, type, numComponents, comp);
403            } else {
404                _dataSet->getDataRange(r, name, type, comp);
405            }
406            _lut->SetRange(r);
407        } else {
408            switch (mode) {
409            case COLOR_CONSTANT:
410            case COLOR_BY_SCALAR:
411                _lut->SetRange(_dataRange);
412                break;
413            case COLOR_BY_VECTOR_MAGNITUDE:
414                _lut->SetRange(_vectorMagnitudeRange);
415                break;
416            case COLOR_BY_VECTOR_X:
417                _lut->SetRange(_vectorComponentRange[0]);
418                break;
419            case COLOR_BY_VECTOR_Y:
420                _lut->SetRange(_vectorComponentRange[1]);
421                break;
422            case COLOR_BY_VECTOR_Z:
423                _lut->SetRange(_vectorComponentRange[2]);
424                break;
425            default:
426                break;
427            }
428        }
429    }
430
431    switch (mode) {
432    case COLOR_BY_SCALAR:
433        _mapper->ScalarVisibilityOn();
434        break;
435    case COLOR_BY_VECTOR_MAGNITUDE:
436        _mapper->ScalarVisibilityOn();
437        if (_lut != NULL) {
438            _lut->SetVectorModeToMagnitude();
439        }
440        break;
441    case COLOR_BY_VECTOR_X:
442        _mapper->ScalarVisibilityOn();
443        if (_lut != NULL) {
444            _lut->SetVectorModeToComponent();
445            _lut->SetVectorComponent(0);
446        }
447        break;
448    case COLOR_BY_VECTOR_Y:
449        _mapper->ScalarVisibilityOn();
450        if (_lut != NULL) {
451            _lut->SetVectorModeToComponent();
452            _lut->SetVectorComponent(1);
453        }
454        break;
455    case COLOR_BY_VECTOR_Z:
456        _mapper->ScalarVisibilityOn();
457        if (_lut != NULL) {
458            _lut->SetVectorModeToComponent();
459            _lut->SetVectorComponent(2);
460        }
461        break;
462    case COLOR_CONSTANT:
463    default:
464        _mapper->ScalarVisibilityOff();
465        break;
466    }
467}
468
469/**
470 * \brief Called when the color map has been edited
471 */
472void PseudoColor::updateColorMap()
473{
474    setColorMap(_colorMap);
475}
476
477/**
478 * \brief Associate a colormap lookup table with the DataSet
479 */
480void PseudoColor::setColorMap(ColorMap *cmap)
481{
482    if (cmap == NULL)
483        return;
484
485    _colorMap = cmap;
486 
487    if (_lut == NULL) {
488        _lut = vtkSmartPointer<vtkLookupTable>::New();
489        if (_mapper != NULL) {
490            _mapper->UseLookupTableScalarRangeOn();
491            _mapper->SetLookupTable(_lut);
492        }
493        _lut->DeepCopy(cmap->getLookupTable());
494        switch (_colorMode) {
495        case COLOR_CONSTANT:
496        case COLOR_BY_SCALAR:
497            _lut->SetRange(_dataRange);
498            break;
499        case COLOR_BY_VECTOR_MAGNITUDE:
500            _lut->SetRange(_vectorMagnitudeRange);
501            break;
502        case COLOR_BY_VECTOR_X:
503            _lut->SetRange(_vectorComponentRange[0]);
504            break;
505        case COLOR_BY_VECTOR_Y:
506            _lut->SetRange(_vectorComponentRange[1]);
507            break;
508        case COLOR_BY_VECTOR_Z:
509            _lut->SetRange(_vectorComponentRange[2]);
510            break;
511        default:
512            break;
513        }
514    } else {
515        double range[2];
516        _lut->GetTableRange(range);
517        _lut->DeepCopy(cmap->getLookupTable());
518        _lut->SetRange(range);
519        _lut->Modified();
520    }
521
522    switch (_colorMode) {
523    case COLOR_BY_VECTOR_MAGNITUDE:
524        _lut->SetVectorModeToMagnitude();
525        break;
526    case COLOR_BY_VECTOR_X:
527        _lut->SetVectorModeToComponent();
528        _lut->SetVectorComponent(0);
529        break;
530    case COLOR_BY_VECTOR_Y:
531        _lut->SetVectorModeToComponent();
532        _lut->SetVectorComponent(1);
533        break;
534    case COLOR_BY_VECTOR_Z:
535        _lut->SetVectorModeToComponent();
536        _lut->SetVectorComponent(2);
537        break;
538    default:
539         break;
540    }
541}
542
543/**
544 * \brief Set a group of world coordinate planes to clip rendering
545 *
546 * Passing NULL for planes will remove all cliping planes
547 */
548void PseudoColor::setClippingPlanes(vtkPlaneCollection *planes)
549{
550    if (_mapper != NULL) {
551        _mapper->SetClippingPlanes(planes);
552    }
553}
Note: See TracBrowser for help on using the repository browser.