source: vtkvis/tags/1.7.2/PseudoColor.cpp @ 4832

Last change on this file since 4832 was 3818, checked in by ldelgass, 11 years ago

Add colormode for PolyData?, by default color map if an active scalar field with
multiple components is present (color scalars). Add protocol for image/text3d
drawing components. Make 'color' synonym for 'ccolor' in protocol commands,
also accept 'ccolor' or 'constant' as colormode. Add experimental commands for
simple legend (no field info required), extended colormode for heightmap with
explicit range. Add fontconfig support (needs testing). Added code to remove
duplicate points from unstructured grids. Should probably be optional. Stack
trace on error (if WANT_TRACE is on) using VTK.

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