source: vtkvis/trunk/HeightMap.cpp @ 4643

Last change on this file since 4643 was 4037, checked in by ldelgass, 10 years ago

Fix square aspect scaling problem

  • Property svn:eol-style set to native
File size: 43.6 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 <vtkDataSet.h>
11#include <vtkPointData.h>
12#include <vtkCellData.h>
13#include <vtkCellDataToPointData.h>
14#include <vtkPolyDataMapper.h>
15#include <vtkUnstructuredGrid.h>
16#include <vtkProperty.h>
17#include <vtkImageData.h>
18#include <vtkLookupTable.h>
19#include <vtkTransform.h>
20#include <vtkDelaunay2D.h>
21#include <vtkDelaunay3D.h>
22#include <vtkGaussianSplatter.h>
23#include <vtkExtractVOI.h>
24#include <vtkDataSetSurfaceFilter.h>
25#include <vtkContourFilter.h>
26#include <vtkStripper.h>
27#include <vtkWarpScalar.h>
28#include <vtkPropAssembly.h>
29#include <vtkCutter.h>
30#include <vtkPlane.h>
31
32#include "HeightMap.h"
33#include "Renderer.h"
34#include "Trace.h"
35
36using namespace VtkVis;
37
38#define TRANSLATE_TO_ZORIGIN
39
40HeightMap::HeightMap(int numContours, double heightScale) :
41    GraphicsObject(),
42    _numContours(numContours),
43    _contourColorMap(false),
44    _contourEdgeWidth(1.0),
45    _warpScale(heightScale),
46    _sliceAxis(Z_AXIS),
47    _pipelineInitialized(false),
48    _cloudStyle(CLOUD_MESH),
49    _colorMap(NULL),
50    _colorMode(COLOR_BY_SCALAR),
51    _colorFieldType(DataSet::POINT_DATA),
52    _renderer(NULL)
53{
54    _contourEdgeColor[0] = 0.0f;
55    _contourEdgeColor[1] = 0.0f;
56    _contourEdgeColor[2] = 0.0f;
57    _colorFieldRange[0] = DBL_MAX;
58    _colorFieldRange[1] = -DBL_MAX;
59}
60
61HeightMap::HeightMap(const std::vector<double>& contours, double heightScale) :
62    GraphicsObject(),
63    _numContours(contours.size()),
64    _contours(contours),
65    _contourColorMap(false),
66    _contourEdgeWidth(1.0),
67    _warpScale(heightScale),
68    _sliceAxis(Z_AXIS),
69    _pipelineInitialized(false),
70    _cloudStyle(CLOUD_MESH),
71    _colorMap(NULL),
72    _colorMode(COLOR_BY_SCALAR),
73    _colorFieldType(DataSet::POINT_DATA),
74    _renderer(NULL)
75{
76    _contourEdgeColor[0] = 0.0f;
77    _contourEdgeColor[1] = 0.0f;
78    _contourEdgeColor[2] = 0.0f;
79    _colorFieldRange[0] = DBL_MAX;
80    _colorFieldRange[1] = -DBL_MAX;
81}
82
83HeightMap::~HeightMap()
84{
85#ifdef WANT_TRACE
86    if (_dataSet != NULL)
87        TRACE("Deleting HeightMap for %s", _dataSet->getName().c_str());
88    else
89        TRACE("Deleting HeightMap with NULL DataSet");
90#endif
91}
92
93void HeightMap::setDataSet(DataSet *dataSet,
94                           Renderer *renderer)
95{
96    if (_dataSet != dataSet) {
97        _dataSet = dataSet;
98        _renderer = renderer;
99
100        if (_dataSet != NULL) {
101            TRACE("DataSet name: '%s' type: %s",
102                  _dataSet->getName().c_str(),
103                  _dataSet->getVtkType());
104
105            if (renderer->getUseCumulativeRange()) {
106                renderer->getCumulativeDataRange(_dataRange,
107                                                 _dataSet->getActiveScalarsName(),
108                                                 1);
109                const char *activeVectors = _dataSet->getActiveVectorsName();
110                if (activeVectors != NULL) {
111                    renderer->getCumulativeDataRange(_vectorMagnitudeRange,
112                                                     activeVectors,
113                                                     3);
114                    for (int i = 0; i < 3; i++) {
115                        renderer->getCumulativeDataRange(_vectorComponentRange[i],
116                                                         activeVectors,
117                                                         3, i);
118                    }
119                }
120            } else {
121                _dataSet->getScalarRange(_dataRange);
122                _dataSet->getVectorRange(_vectorMagnitudeRange);
123                for (int i = 0; i < 3; i++) {
124                    _dataSet->getVectorRange(_vectorComponentRange[i], i);
125                }
126            }
127        }
128
129        update();
130    }
131}
132
133/**
134 * Compute a data scaling factor to make maximum height (at default _warpScale)
135 * equivalent to largest dimension of data set
136 */
137void HeightMap::computeDataScale()
138{
139    if (_dataSet == NULL)
140        return;
141
142    double bounds[6];
143    double boundsRange = 0.0;
144    _dataSet->getBounds(bounds);
145    for (int i = 0; i < 6; i+=2) {
146        double r = bounds[i+1] - bounds[i];
147        if (r > boundsRange)
148            boundsRange = r;
149    }
150    double datRange = _dataRange[1] - _dataRange[0];
151    if (datRange < 1.0e-17) {
152        _dataScale = 1.0;
153    } else {
154        _dataScale = boundsRange / datRange;
155    }
156    TRACE("Bounds range: %g data range: %g scaling: %g",
157          boundsRange, datRange, _dataScale);
158}
159
160/**
161 * \brief Create and initialize VTK Props to render the colormapped dataset
162 */
163void HeightMap::initProp()
164{
165    if (_dsActor == NULL) {
166        _dsActor = vtkSmartPointer<vtkActor>::New();
167        _dsActor->GetProperty()->SetOpacity(_opacity);
168        _dsActor->GetProperty()->SetColor(_color[0],
169                                          _color[1],
170                                          _color[2]);
171        _dsActor->GetProperty()->SetEdgeColor(_edgeColor[0],
172                                              _edgeColor[1],
173                                              _edgeColor[2]);
174        _dsActor->GetProperty()->SetLineWidth(_edgeWidth);
175        _dsActor->GetProperty()->EdgeVisibilityOff();
176        _dsActor->GetProperty()->SetAmbient(0.2);
177        _dsActor->GetProperty()->LightingOn();
178    }
179    if (_contourActor == NULL) {
180        _contourActor = vtkSmartPointer<vtkActor>::New();
181        _contourActor->GetProperty()->SetOpacity(_opacity);
182        _contourActor->GetProperty()->SetColor(_contourEdgeColor[0],
183                                               _contourEdgeColor[1],
184                                               _contourEdgeColor[2]);
185        _contourActor->GetProperty()->SetEdgeColor(_contourEdgeColor[0],
186                                                   _contourEdgeColor[1],
187                                                   _contourEdgeColor[2]);
188        _contourActor->GetProperty()->SetLineWidth(_contourEdgeWidth);
189        _contourActor->GetProperty()->EdgeVisibilityOff();
190        _contourActor->GetProperty()->SetAmbient(1.0);
191        _contourActor->GetProperty()->LightingOff();
192    }
193    // Contour actor is added in update() if contours are produced
194    if (_prop == NULL) {
195        _prop = vtkSmartPointer<vtkAssembly>::New();
196        getAssembly()->AddPart(_dsActor);
197    } else {
198        getAssembly()->RemovePart(_contourActor);
199    }
200}
201
202/**
203 * \brief Internal method to set up pipeline after a state change
204 */
205void HeightMap::update()
206{
207    if (_dataSet == NULL)
208        return;
209
210    TRACE("DataSet: %s", _dataSet->getName().c_str());
211
212    vtkDataSet *ds = _dataSet->getVtkDataSet();
213
214    computeDataScale();
215
216    // Contour filter to generate isolines
217    if (_contourFilter == NULL) {
218        _contourFilter = vtkSmartPointer<vtkContourFilter>::New();
219    }
220
221    // Mapper, actor to render color-mapped data set
222    if (_mapper == NULL) {
223        _mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
224        // Map scalars through lookup table regardless of type
225        _mapper->SetColorModeToMapScalars();
226    }
227
228    if (!_pipelineInitialized) {
229        vtkSmartPointer<vtkCellDataToPointData> cellToPtData;
230
231        if (ds->GetPointData() == NULL ||
232            ds->GetPointData()->GetScalars() == NULL) {
233            TRACE("No scalar point data in dataset %s", _dataSet->getName().c_str());
234            if (ds->GetCellData() != NULL &&
235                ds->GetCellData()->GetScalars() != NULL) {
236                cellToPtData =
237                    vtkSmartPointer<vtkCellDataToPointData>::New();
238#ifdef USE_VTK6
239                cellToPtData->SetInputData(ds);
240#else
241                cellToPtData->SetInput(ds);
242#endif
243                //cellToPtData->PassCellDataOn();
244                cellToPtData->Update();
245                ds = cellToPtData->GetOutput();
246            } else {
247                USER_ERROR("No scalar field was found in the data set.");
248                return;
249            }
250        }
251
252        _splatter = NULL;
253
254        if (_dataSet->isCloud()) {
255            // DataSet is a point cloud
256            PrincipalPlane plane;
257            double offset;
258            if (_dataSet->is2D(&plane, &offset)) {
259                if (_cloudStyle == CLOUD_MESH) {
260                    // Result of Delaunay2D is a PolyData
261                    vtkSmartPointer<vtkDelaunay2D> mesher = vtkSmartPointer<vtkDelaunay2D>::New();
262                    if (plane == PLANE_ZY) {
263                        vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
264                        trans->RotateWXYZ(90, 0, 1, 0);
265                        if (offset != 0.0) {
266                            trans->Translate(-offset, 0, 0);
267                        }
268                        mesher->SetTransform(trans);
269                        _sliceAxis = X_AXIS;
270                    } else if (plane == PLANE_XZ) {
271                        vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
272                        trans->RotateWXYZ(-90, 1, 0, 0);
273                        if (offset != 0.0) {
274                            trans->Translate(0, -offset, 0);
275                        }
276                        mesher->SetTransform(trans);
277                        _sliceAxis = Y_AXIS;
278                    } else if (offset != 0.0) {
279                        // XY with Z offset
280                        vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
281                        trans->Translate(0, 0, -offset);
282                        mesher->SetTransform(trans);
283                    }
284#ifdef USE_VTK6
285                    mesher->SetInputData(ds);
286#else
287                    mesher->SetInput(ds);
288#endif
289                    vtkAlgorithmOutput *warpOutput = initWarp(mesher->GetOutputPort());
290                    _mapper->SetInputConnection(warpOutput);
291                    _contourFilter->SetInputConnection(warpOutput);
292                } else {
293                    // _cloudStyle == CLOUD_SPLAT
294                    if (_splatter == NULL)
295                        _splatter = vtkSmartPointer<vtkGaussianSplatter>::New();
296                    if (_volumeSlicer == NULL)
297                        _volumeSlicer = vtkSmartPointer<vtkExtractVOI>::New();
298#ifdef USE_VTK6
299                    _splatter->SetInputData(ds);
300#else
301                    _splatter->SetInput(ds);
302#endif
303                    int dims[3];
304                    _splatter->GetSampleDimensions(dims);
305                    TRACE("Sample dims: %d %d %d", dims[0], dims[1], dims[2]);
306                    if (plane == PLANE_ZY) {
307                        dims[0] = 3;
308                        _volumeSlicer->SetVOI(1, 1, 0, dims[1]-1, 0, dims[1]-1);
309                        _sliceAxis = X_AXIS;
310                    } else if (plane == PLANE_XZ) {
311                        dims[1] = 3;
312                        _volumeSlicer->SetVOI(0, dims[0]-1, 1, 1, 0, dims[2]-1);
313                        _sliceAxis = Y_AXIS;
314                    } else {
315                        dims[2] = 3;
316                        _volumeSlicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, 1, 1);
317                    }
318                    _splatter->SetSampleDimensions(dims);
319                    double bounds[6];
320                    _splatter->Update();
321                    _splatter->GetModelBounds(bounds);
322                    TRACE("Model bounds: %g %g %g %g %g %g",
323                          bounds[0], bounds[1],
324                          bounds[2], bounds[3],
325                          bounds[4], bounds[5]);
326                    _volumeSlicer->SetInputConnection(_splatter->GetOutputPort());
327                    _volumeSlicer->SetSampleRate(1, 1, 1);
328                    vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
329                    gf->UseStripsOn();
330                    gf->SetInputConnection(_volumeSlicer->GetOutputPort());
331                    vtkAlgorithmOutput *warpOutput = initWarp(gf->GetOutputPort());
332                    _mapper->SetInputConnection(warpOutput);
333                    _contourFilter->SetInputConnection(warpOutput);
334                }
335            } else {
336                // 3D point cloud
337                vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
338                gf->UseStripsOn();
339                if (_cloudStyle == CLOUD_MESH) {
340                     // Result of Delaunay3D mesher is unstructured grid
341                    vtkSmartPointer<vtkDelaunay3D> mesher = vtkSmartPointer<vtkDelaunay3D>::New();
342#ifdef USE_VTK6
343                    mesher->SetInputData(ds);
344#else
345                    mesher->SetInput(ds);
346#endif
347                    // Run the mesher
348                    mesher->Update();
349                    // Get bounds of resulting grid
350                    double bounds[6];
351                    mesher->GetOutput()->GetBounds(bounds);
352                    // Sample a plane within the grid bounding box
353                    vtkSmartPointer<vtkCutter> cutter = vtkSmartPointer<vtkCutter>::New();
354                    cutter->SetInputConnection(mesher->GetOutputPort());
355                    if (_cutPlane == NULL) {
356                        _cutPlane = vtkSmartPointer<vtkPlane>::New();
357                    }
358                    _cutPlane->SetNormal(0, 0, 1);
359                    _cutPlane->SetOrigin(0,
360                                         0,
361                                         bounds[4] + (bounds[5]-bounds[4])/2.);
362                    cutter->SetCutFunction(_cutPlane);
363                    gf->SetInputConnection(cutter->GetOutputPort());
364                } else {
365                    // _cloudStyle == CLOUD_SPLAT
366                    if (_splatter == NULL)
367                        _splatter = vtkSmartPointer<vtkGaussianSplatter>::New();
368#ifdef USE_VTK6
369                    _splatter->SetInputData(ds);
370#else
371                    _splatter->SetInput(ds);
372#endif
373                    int dims[3];
374                    _splatter->GetSampleDimensions(dims);
375                    TRACE("Sample dims: %d %d %d", dims[0], dims[1], dims[2]);
376                    dims[2] = 3;
377                    _splatter->SetSampleDimensions(dims);
378                    double bounds[6];
379                    _splatter->Update();
380                    _splatter->GetModelBounds(bounds);
381                    TRACE("Model bounds: %g %g %g %g %g %g",
382                          bounds[0], bounds[1],
383                          bounds[2], bounds[3],
384                          bounds[4], bounds[5]);
385                    if (_volumeSlicer == NULL)
386                        _volumeSlicer = vtkSmartPointer<vtkExtractVOI>::New();
387                    _volumeSlicer->SetInputConnection(_splatter->GetOutputPort());
388                    _volumeSlicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, 1, 1);
389                    _volumeSlicer->SetSampleRate(1, 1, 1);
390                    gf->SetInputConnection(_volumeSlicer->GetOutputPort());
391                }
392                vtkAlgorithmOutput *warpOutput = initWarp(gf->GetOutputPort());
393                _mapper->SetInputConnection(warpOutput);
394                _contourFilter->SetInputConnection(warpOutput);
395            }
396        } else if (vtkPolyData::SafeDownCast(ds) != NULL) {
397            // DataSet is a vtkPolyData with lines and/or polygons
398            vtkAlgorithmOutput *warpOutput = initWarp(ds);
399            if (warpOutput != NULL) {
400                _mapper->SetInputConnection(warpOutput);
401                _contourFilter->SetInputConnection(warpOutput);
402            } else {
403#ifdef USE_VTK6
404                _mapper->SetInputData(vtkPolyData::SafeDownCast(ds));
405                _contourFilter->SetInputData(ds);
406#else
407                _mapper->SetInput(vtkPolyData::SafeDownCast(ds));
408                _contourFilter->SetInput(ds);
409#endif
410            }
411        } else {
412            // DataSet is NOT a vtkPolyData
413            // Can be: image/volume/uniform grid, structured grid, unstructured grid, rectilinear grid
414            vtkSmartPointer<vtkDataSetSurfaceFilter> gf = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
415            gf->UseStripsOn();
416            vtkImageData *imageData = vtkImageData::SafeDownCast(ds);
417            if (!_dataSet->is2D() && imageData != NULL) {
418                // 3D image/volume/uniform grid
419                if (_volumeSlicer == NULL)
420                    _volumeSlicer = vtkSmartPointer<vtkExtractVOI>::New();
421                int dims[3];
422                imageData->GetDimensions(dims);
423                TRACE("Image data dimensions: %d %d %d", dims[0], dims[1], dims[2]);
424#ifdef USE_VTK6
425                _volumeSlicer->SetInputData(ds);
426#else
427                _volumeSlicer->SetInput(ds);
428#endif
429                _volumeSlicer->SetVOI(0, dims[0]-1, 0, dims[1]-1, (dims[2]-1)/2, (dims[2]-1)/2);
430                _volumeSlicer->SetSampleRate(1, 1, 1);
431                gf->SetInputConnection(_volumeSlicer->GetOutputPort());
432            } else if (!_dataSet->is2D() && imageData == NULL) {
433                // 3D structured grid, unstructured grid, or rectilinear grid
434                double bounds[6];
435                ds->GetBounds(bounds);
436                // Sample a plane within the grid bounding box
437                vtkSmartPointer<vtkCutter> cutter = vtkSmartPointer<vtkCutter>::New();
438#ifdef USE_VTK6
439                cutter->SetInputData(ds);
440#else
441                cutter->SetInput(ds);
442#endif
443                if (_cutPlane == NULL) {
444                    _cutPlane = vtkSmartPointer<vtkPlane>::New();
445                }
446                _cutPlane->SetNormal(0, 0, 1);
447                _cutPlane->SetOrigin(0,
448                                     0,
449                                     bounds[4] + (bounds[5]-bounds[4])/2.);
450                cutter->SetCutFunction(_cutPlane);
451                gf->SetInputConnection(cutter->GetOutputPort());
452            } else {
453                // 2D data
454#ifdef USE_VTK6
455                gf->SetInputData(ds);
456#else
457                gf->SetInput(ds);
458#endif
459            }
460            vtkAlgorithmOutput *warpOutput = initWarp(gf->GetOutputPort());
461            _mapper->SetInputConnection(warpOutput);
462            _contourFilter->SetInputConnection(warpOutput);
463        }
464    }
465
466    _contourFilter->ComputeNormalsOff();
467    _contourFilter->ComputeGradientsOff();
468
469    // Speed up multiple contour computation at cost of extra memory use
470    if (_numContours > 1) {
471        _contourFilter->UseScalarTreeOn();
472    } else {
473        _contourFilter->UseScalarTreeOff();
474    }
475
476    _contourFilter->SetNumberOfContours(_numContours);
477
478    if (_contours.empty()) {
479        // Evenly spaced isovalues
480        if (_colorFieldRange[0] <= _colorFieldRange[1]) {
481            _contourFilter->GenerateValues(_numContours, _colorFieldRange[0], _colorFieldRange[1]);
482        } else {
483            _contourFilter->GenerateValues(_numContours, _dataRange[0], _dataRange[1]);
484        }
485    } else {
486        // User-supplied isovalues
487        for (int i = 0; i < _numContours; i++) {
488            _contourFilter->SetValue(i, _contours[i]);
489        }
490    }
491
492    initProp();
493
494    if (_contourMapper == NULL) {
495        _contourMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
496        _contourMapper->ScalarVisibilityOff();
497        _contourMapper->SetResolveCoincidentTopologyToPolygonOffset();
498        vtkSmartPointer<vtkStripper> stripper = vtkSmartPointer<vtkStripper>::New();
499        stripper->SetInputConnection(_contourFilter->GetOutputPort());
500        _contourMapper->SetInputConnection(stripper->GetOutputPort());
501        _contourActor->SetMapper(_contourMapper);
502        //_contourActor->InterpolateScalarsBeforeMappingOn();
503    }
504
505    setInterpolateBeforeMapping(true);
506
507    if (_lut == NULL) {
508        setColorMap(ColorMap::getDefault());
509        setColorMode(_colorMode);
510    } else if (!_pipelineInitialized) {
511        double *rangePtr = _colorFieldRange;
512        if (_colorFieldRange[0] > _colorFieldRange[1]) {
513            rangePtr = NULL;
514        }
515        setColorMode(_colorMode, _colorFieldType, _colorFieldName.c_str(), rangePtr);
516    }
517
518    _pipelineInitialized = true;
519
520    // Ensure updated dataScale is applied
521    if (_warp != NULL) {
522        _warp->SetScaleFactor(_warpScale * _dataScale);
523    }
524#ifdef TRANSLATE_TO_ZORIGIN
525    double pos[3];
526    pos[0] = 0;
527    pos[1] = 0;
528    pos[2] = - _warpScale * _dataScale * _dataRange[0];
529    setPosition(pos);
530    TRACE("Z_POS: %g", pos[2]);
531#endif
532    _dsActor->SetMapper(_mapper);
533
534    _mapper->Update();
535    _contourMapper->Update();
536
537    // Only add contour actor to assembly if contours were
538    // produced, in order to prevent messing up assembly bounds
539    double bounds[6];
540    _contourActor->GetBounds(bounds);
541    if (bounds[0] <= bounds[1]) {
542        getAssembly()->AddPart(_contourActor);
543    }
544}
545
546void HeightMap::setCloudStyle(CloudStyle style)
547{
548    if (style != _cloudStyle) {
549        _cloudStyle = style;
550        if (_dataSet != NULL) {
551            _pipelineInitialized = false;
552            update();
553        }
554    }
555}
556
557void HeightMap::setInterpolateBeforeMapping(bool state)
558{
559    if (_mapper != NULL) {
560        _mapper->SetInterpolateScalarsBeforeMapping((state ? 1 : 0));
561    }
562}
563
564void HeightMap::setAspect(double aspect)
565{
566    double bounds[6];
567    vtkDataSet *ds = _dataSet->getVtkDataSet();
568    ds->GetBounds(bounds);
569    double size[3];
570    size[0] = bounds[1] - bounds[0];
571    size[1] = bounds[3] - bounds[2];
572    size[2] = bounds[5] - bounds[4];
573    double scale[3];
574    scale[0] = scale[1] = scale[2] = 1.;
575
576    if (aspect == 1.0) {
577        // Square
578        switch (_sliceAxis) {
579        case X_AXIS: {
580            if (size[1] > size[2] && size[2] > 0.0) {
581                scale[2] = size[1] / size[2];
582            } else if (size[2] > size[1] && size[1] > 0.0) {
583                scale[1] = size[2] / size[1];
584            }
585        }
586            break;
587        case Y_AXIS: {
588            if (size[0] > size[2] && size[2] > 0.0) {
589                scale[2] = size[0] / size[2];
590            } else if (size[2] > size[0] && size[0] > 0.0) {
591                scale[0] = size[2] / size[0];
592            }
593        }
594            break;
595        case Z_AXIS: {
596            if (size[0] > size[1] && size[1] > 0.0) {
597                scale[1] = size[0] / size[1];
598            } else if (size[1] > size[0] && size[0] > 0.0) {
599                scale[0] = size[1] / size[0];
600            }
601        }
602            break;
603        }
604    } else if (aspect != 0.0) {
605        switch (_sliceAxis) {
606        case X_AXIS: {
607            if (aspect > 1.0) {
608                if (size[2] > size[1] && size[1] > 0.0) {
609                    scale[1] = (size[2] / aspect) / size[1];
610                } else if (size[2] > 0.0) {
611                    scale[2] = (size[1] * aspect) / size[2];
612                }
613            } else {
614                if (size[1] > size[2] && size[2] > 0.0) {
615                    scale[2] = (size[1] * aspect) / size[2];
616                } else if (size[1] > 0.0) {
617                    scale[1] = (size[2] / aspect) / size[1];
618                }
619            }
620        }
621            break;
622        case Y_AXIS: {
623            if (aspect > 1.0) {
624                if (size[0] > size[2] && size[2] > 0.0) {
625                    scale[2] = (size[0] / aspect) / size[2];
626                } else if (size[0] > 0.0) {
627                    scale[0] = (size[2] * aspect) / size[0];
628                }
629            } else {
630                if (size[2] > size[0] && size[0] > 0.0) {
631                    scale[0] = (size[2] * aspect) / size[0];
632                } else if (size[2] > 0.0) {
633                    scale[2] = (size[0] / aspect) / size[2];
634                }
635            }
636        }
637            break;
638        case Z_AXIS: {
639            if (aspect > 1.0) {
640                if (size[0] > size[1] && size[1] > 0.0) {
641                    scale[1] = (size[0] / aspect) / size[1];
642                } else if (size[0] > 0.0) {
643                    scale[0] = (size[1] * aspect) / size[0];
644                }
645            } else {
646                if (size[1] > size[0] && size[0] > 0.0) {
647                    scale[0] = (size[1] * aspect) / size[0];
648                } else if (size[1] > 0.0) {
649                    scale[1] = (size[0] / aspect) / size[1];
650                }
651            }
652        }
653            break;
654        }
655    }
656
657    TRACE("%s dims %g,%g,%g", _dataSet->getName().c_str(),
658          size[0], size[1], size[2]);
659    TRACE("Setting scale to %g,%g,%g", scale[0], scale[1], scale[2]);
660    setScale(scale);
661}
662
663vtkAlgorithmOutput *HeightMap::initWarp(vtkAlgorithmOutput *input)
664{
665    TRACE("Warp scale: %g", _warpScale);
666    if (_warpScale == 0.0) {
667        _warp = NULL;
668        _normalsGenerator = NULL;
669        return input;
670    } else if (input == NULL) {
671        ERROR("NULL input");
672        return input;
673    } else {
674        if (_warp == NULL)
675            _warp = vtkSmartPointer<vtkWarpScalar>::New();
676        if (_normalsGenerator == NULL)
677            _normalsGenerator = vtkSmartPointer<vtkPolyDataNormals>::New();
678        switch (_sliceAxis) {
679        case X_AXIS:
680            _warp->SetNormal(1, 0, 0);
681            break;
682        case Y_AXIS:
683            _warp->SetNormal(0, 1, 0);
684            break;
685        default:
686            _warp->SetNormal(0, 0, 1);
687        }
688        _warp->UseNormalOn();
689        _warp->SetScaleFactor(_warpScale * _dataScale);
690        _warp->SetInputConnection(input);
691        _normalsGenerator->SetInputConnection(_warp->GetOutputPort());
692        _normalsGenerator->SetFeatureAngle(90.);
693        _normalsGenerator->AutoOrientNormalsOff();
694        _normalsGenerator->ComputePointNormalsOn();
695        return _normalsGenerator->GetOutputPort();
696    }
697}
698
699vtkAlgorithmOutput *HeightMap::initWarp(vtkDataSet *dsInput)
700{
701    TRACE("Warp scale: %g", _warpScale);
702    if (_warpScale == 0.0) {
703        _warp = NULL;
704        _normalsGenerator = NULL;
705        return NULL;
706    } else if (dsInput == NULL) {
707        ERROR("NULL input");
708        return NULL;
709    } else {
710        if (_warp == NULL)
711            _warp = vtkSmartPointer<vtkWarpScalar>::New();
712        if (_normalsGenerator == NULL)
713            _normalsGenerator = vtkSmartPointer<vtkPolyDataNormals>::New();
714        switch (_sliceAxis) {
715        case X_AXIS:
716            _warp->SetNormal(1, 0, 0);
717            break;
718        case Y_AXIS:
719            _warp->SetNormal(0, 1, 0);
720            break;
721        default:
722            _warp->SetNormal(0, 0, 1);
723        }
724        _warp->UseNormalOn();
725        _warp->SetScaleFactor(_warpScale * _dataScale);
726#ifdef USE_VTK6
727        _warp->SetInputData(dsInput);
728#else
729        _warp->SetInput(dsInput);
730#endif
731        _normalsGenerator->SetInputConnection(_warp->GetOutputPort());
732        _normalsGenerator->SetFeatureAngle(90.);
733        _normalsGenerator->AutoOrientNormalsOff();
734        _normalsGenerator->ComputePointNormalsOn();
735        return _normalsGenerator->GetOutputPort();
736    }
737}
738
739/**
740 * \brief Controls relative scaling of height of mountain plot
741 */
742void HeightMap::setHeightScale(double scale)
743{
744    if (_warpScale == scale)
745        return;
746
747    _warpScale = scale;
748    if (_warp == NULL) {
749        vtkAlgorithmOutput *warpOutput = initWarp(_mapper->GetInputConnection(0, 0));
750        _mapper->SetInputConnection(warpOutput);
751        _contourFilter->SetInputConnection(warpOutput);
752    } else if (scale == 0.0) {
753        vtkAlgorithmOutput *warpInput = _warp->GetInputConnection(0, 0);
754        _mapper->SetInputConnection(warpInput);
755        _contourFilter->SetInputConnection(warpInput);
756        _warp = NULL;
757    } else {
758        _warp->SetScaleFactor(_warpScale * _dataScale);
759    }
760#ifdef TRANSLATE_TO_ZORIGIN
761    double pos[3];
762    pos[0] = 0;
763    pos[1] = 0;
764    pos[2] = - _warpScale * _dataScale * _dataRange[0];
765    setPosition(pos);
766    TRACE("Z_POS: %g", pos[2]);
767#endif
768    if (_mapper != NULL)
769        _mapper->Update();
770    if (_contourMapper != NULL)
771        _contourMapper->Update();
772}
773
774/**
775 * \brief Select a 2D slice plane from a 3D DataSet
776 *
777 * \param[in] axis Axis of slice plane
778 * \param[in] ratio Position [0,1] of slice plane along axis
779 */
780void HeightMap::selectVolumeSlice(Axis axis, double ratio)
781{
782    if (_dataSet->is2D()) {
783        WARN("DataSet not 3D, returning");
784        return;
785    }
786
787    if (_volumeSlicer == NULL &&
788        _cutPlane == NULL) {
789        WARN("Called before update() or DataSet is not a volume");
790        return;
791    }
792
793    _sliceAxis = axis;
794    if (_warp != NULL) {
795        switch (axis) {
796        case X_AXIS:
797            _warp->SetNormal(1, 0, 0);
798            break;
799        case Y_AXIS:
800            _warp->SetNormal(0, 1, 0);
801            break;
802        case Z_AXIS:
803            _warp->SetNormal(0, 0, 1);
804            break;
805        default:
806            ERROR("Invalid Axis");
807            return;
808        }
809    }
810
811    if (_cutPlane != NULL) {
812        double bounds[6];
813        _dataSet->getBounds(bounds);
814        switch (axis) {
815        case X_AXIS:
816            _cutPlane->SetNormal(1, 0, 0);
817            _cutPlane->SetOrigin(bounds[0] + (bounds[1]-bounds[0]) * ratio,
818                                 0,
819                                 0);
820            break;
821        case Y_AXIS:
822            _cutPlane->SetNormal(0, 1, 0);
823            _cutPlane->SetOrigin(0,
824                                 bounds[2] + (bounds[3]-bounds[2]) * ratio,
825                                 0);
826            break;
827        case Z_AXIS:
828            _cutPlane->SetNormal(0, 0, 1);
829            _cutPlane->SetOrigin(0,
830                                 0,
831                                 bounds[4] + (bounds[5]-bounds[4]) * ratio);
832            break;
833        default:
834            ERROR("Invalid Axis");
835            return;
836        }
837    } else {
838        int dims[3];
839        if (_splatter != NULL) {
840            _splatter->GetSampleDimensions(dims);
841        } else {
842            vtkImageData *imageData = vtkImageData::SafeDownCast(_dataSet->getVtkDataSet());
843            if (imageData == NULL) {
844                ERROR("Not a volume data set");
845                return;
846            }
847            imageData->GetDimensions(dims);
848        }
849        int voi[6];
850
851        switch (axis) {
852        case X_AXIS:
853            voi[0] = voi[1] = (int)((dims[0]-1) * ratio);
854            voi[2] = 0;
855            voi[3] = dims[1]-1;
856            voi[4] = 0;
857            voi[5] = dims[2]-1;
858            break;
859        case Y_AXIS:
860            voi[0] = 0;
861            voi[1] = dims[0]-1;
862            voi[2] = voi[3] = (int)((dims[1]-1) * ratio);
863            voi[4] = 0;
864            voi[5] = dims[2]-1;
865            break;
866        case Z_AXIS:
867            voi[0] = 0;
868            voi[1] = dims[0]-1;
869            voi[2] = 0;
870            voi[3] = dims[1]-1;
871            voi[4] = voi[5] = (int)((dims[2]-1) * ratio);
872            break;
873        default:
874            ERROR("Invalid Axis");
875            return;
876        }
877
878        _volumeSlicer->SetVOI(voi);
879    }
880
881    if (_mapper != NULL)
882        _mapper->Update();
883    if (_contourMapper != NULL)
884        _contourMapper->Update();
885}
886
887void HeightMap::updateRanges(Renderer *renderer)
888{
889    TRACE("Enter");
890
891    if (_dataSet == NULL) {
892        ERROR("called before setDataSet");
893        return;
894    }
895
896    if (renderer->getUseCumulativeRange()) {
897        renderer->getCumulativeDataRange(_dataRange,
898                                         _dataSet->getActiveScalarsName(),
899                                         1);
900        const char *activeVectors = _dataSet->getActiveVectorsName();
901        if (activeVectors != NULL) {
902            renderer->getCumulativeDataRange(_vectorMagnitudeRange,
903                                             activeVectors,
904                                             3);
905            for (int i = 0; i < 3; i++) {
906                renderer->getCumulativeDataRange(_vectorComponentRange[i],
907                                                 activeVectors,
908                                                 3, i);
909            }
910        }
911    } else {
912        _dataSet->getScalarRange(_dataRange);
913        _dataSet->getVectorRange(_vectorMagnitudeRange);
914        for (int i = 0; i < 3; i++) {
915            _dataSet->getVectorRange(_vectorComponentRange[i], i);
916        }
917    }
918 
919    // Need to update color map ranges
920    double *rangePtr = _colorFieldRange;
921    if (_colorFieldRange[0] > _colorFieldRange[1]) {
922        rangePtr = NULL;
923    }
924    setColorMode(_colorMode, _colorFieldType, _colorFieldName.c_str(), rangePtr);
925
926    if ((_contours.empty() && _numContours > 0) ||
927        _warpScale > 0.0) {
928        // Contour isovalues and/or warp need to be recomputed
929        update();
930    }
931
932    TRACE("Leave");
933}
934
935void HeightMap::setContourLineColorMapEnabled(bool mode)
936{
937    _contourColorMap = mode;
938
939    double *rangePtr = _colorFieldRange;
940    if (_colorFieldRange[0] > _colorFieldRange[1]) {
941        rangePtr = NULL;
942    }
943    setColorMode(_colorMode, _colorFieldType, _colorFieldName.c_str(), rangePtr);
944}
945
946void HeightMap::setColorMode(ColorMode mode)
947{
948    _colorMode = mode;
949    if (_dataSet == NULL)
950        return;
951
952    switch (mode) {
953    case COLOR_BY_SCALAR:
954        setColorMode(mode,
955                     _dataSet->getActiveScalarsType(),
956                     _dataSet->getActiveScalarsName(),
957                     _dataRange);
958        break;
959    case COLOR_BY_VECTOR_MAGNITUDE:
960        setColorMode(mode,
961                     _dataSet->getActiveVectorsType(),
962                     _dataSet->getActiveVectorsName(),
963                     _vectorMagnitudeRange);
964        break;
965    case COLOR_BY_VECTOR_X:
966        setColorMode(mode,
967                     _dataSet->getActiveVectorsType(),
968                     _dataSet->getActiveVectorsName(),
969                     _vectorComponentRange[0]);
970        break;
971    case COLOR_BY_VECTOR_Y:
972        setColorMode(mode,
973                     _dataSet->getActiveVectorsType(),
974                     _dataSet->getActiveVectorsName(),
975                     _vectorComponentRange[1]);
976        break;
977    case COLOR_BY_VECTOR_Z:
978        setColorMode(mode,
979                     _dataSet->getActiveVectorsType(),
980                     _dataSet->getActiveVectorsName(),
981                     _vectorComponentRange[2]);
982        break;
983    case COLOR_CONSTANT:
984    default:
985        setColorMode(mode, DataSet::POINT_DATA, NULL, NULL);
986        break;
987    }
988}
989
990void HeightMap::setColorMode(ColorMode mode,
991                             const char *name, double range[2])
992{
993    if (_dataSet == NULL)
994        return;
995    DataSet::DataAttributeType type = DataSet::POINT_DATA;
996    int numComponents = 1;
997    if (name != NULL && strlen(name) > 0 &&
998        !_dataSet->getFieldInfo(name, &type, &numComponents)) {
999        ERROR("Field not found: %s", name);
1000        return;
1001    }
1002    setColorMode(mode, type, name, range);
1003}
1004
1005void HeightMap::setColorMode(ColorMode mode, DataSet::DataAttributeType type,
1006                             const char *name, double range[2])
1007{
1008    _colorMode = mode;
1009    _colorFieldType = type;
1010    if (name == NULL)
1011        _colorFieldName.clear();
1012    else
1013        _colorFieldName = name;
1014    if (range == NULL) {
1015        _colorFieldRange[0] = DBL_MAX;
1016        _colorFieldRange[1] = -DBL_MAX;
1017    } else {
1018        memcpy(_colorFieldRange, range, sizeof(double)*2);
1019    }
1020
1021    if (_dataSet == NULL || _mapper == NULL)
1022        return;
1023
1024    switch (type) {
1025    case DataSet::POINT_DATA:
1026        _mapper->SetScalarModeToUsePointFieldData();
1027        _contourMapper->SetScalarModeToUsePointFieldData();
1028        break;
1029    case DataSet::CELL_DATA:
1030        _mapper->SetScalarModeToUseCellFieldData();
1031        _contourMapper->SetScalarModeToUseCellFieldData();
1032        break;
1033    default:
1034        ERROR("Unsupported DataAttributeType: %d", type);
1035        return;
1036    }
1037
1038    if (_splatter != NULL) {
1039        if (name != NULL && strlen(name) > 0) {
1040            _splatter->SetInputArrayToProcess(0, 0, 0, vtkDataObject::FIELD_ASSOCIATION_POINTS, name);
1041        }
1042        _mapper->SelectColorArray("SplatterValues");
1043        _contourMapper->SelectColorArray("SplatterValues");
1044    } else if (name != NULL && strlen(name) > 0) {
1045        _mapper->SelectColorArray(name);
1046        _contourMapper->SelectColorArray(name);
1047    } else {
1048        _mapper->SetScalarModeToDefault();
1049        _contourMapper->SetScalarModeToDefault();
1050    }
1051
1052    if (_lut != NULL) {
1053        if (range != NULL) {
1054            _lut->SetRange(range);
1055            TRACE("range %g, %g", range[0], range[1]);
1056        } else if (name != NULL && strlen(name) > 0) {
1057            double r[2];
1058            int comp = -1;
1059            if (mode == COLOR_BY_VECTOR_X)
1060                comp = 0;
1061            else if (mode == COLOR_BY_VECTOR_Y)
1062                comp = 1;
1063            else if (mode == COLOR_BY_VECTOR_Z)
1064                comp = 2;
1065
1066            if (_renderer->getUseCumulativeRange()) {
1067                int numComponents;
1068                if  (!_dataSet->getFieldInfo(name, type, &numComponents)) {
1069                    ERROR("Field not found: %s, type: %d", name, type);
1070                    return;
1071                } else if (numComponents < comp+1) {
1072                    ERROR("Request for component %d in field with %d components",
1073                          comp, numComponents);
1074                    return;
1075                }
1076                _renderer->getCumulativeDataRange(r, name, type, numComponents, comp);
1077            } else {
1078                _dataSet->getDataRange(r, name, type, comp);
1079            }
1080            _lut->SetRange(r);
1081            TRACE("%s range %g, %g", name, r[0], r[1]);
1082        }
1083        _lut->Modified();
1084    }
1085
1086    if (_contourColorMap) {
1087        _contourMapper->ScalarVisibilityOn();
1088    } else {
1089        _contourMapper->ScalarVisibilityOff();
1090    }
1091
1092    switch (mode) {
1093    case COLOR_BY_SCALAR:
1094        _mapper->ScalarVisibilityOn();
1095        break;
1096    case COLOR_BY_VECTOR_MAGNITUDE:
1097        _mapper->ScalarVisibilityOn();
1098        if (_lut != NULL) {
1099            _lut->SetVectorModeToMagnitude();
1100        }
1101        break;
1102    case COLOR_BY_VECTOR_X:
1103        _mapper->ScalarVisibilityOn();
1104        if (_lut != NULL) {
1105            _lut->SetVectorModeToComponent();
1106            _lut->SetVectorComponent(0);
1107        }
1108        break;
1109    case COLOR_BY_VECTOR_Y:
1110        _mapper->ScalarVisibilityOn();
1111        if (_lut != NULL) {
1112            _lut->SetVectorModeToComponent();
1113            _lut->SetVectorComponent(1);
1114        }
1115        break;
1116    case COLOR_BY_VECTOR_Z:
1117        _mapper->ScalarVisibilityOn();
1118        if (_lut != NULL) {
1119            _lut->SetVectorModeToComponent();
1120            _lut->SetVectorComponent(2);
1121        }
1122        break;
1123    case COLOR_CONSTANT:
1124    default:
1125        _mapper->ScalarVisibilityOff();
1126        break;
1127    }
1128}
1129
1130/**
1131 * \brief Called when the color map has been edited
1132 */
1133void HeightMap::updateColorMap()
1134{
1135    setColorMap(_colorMap);
1136}
1137
1138/**
1139 * \brief Associate a colormap lookup table with the DataSet
1140 */
1141void HeightMap::setColorMap(ColorMap *cmap)
1142{
1143    if (cmap == NULL)
1144        return;
1145
1146    _colorMap = cmap;
1147
1148    if (_lut == NULL) {
1149        _lut = vtkSmartPointer<vtkLookupTable>::New();
1150        if (_mapper != NULL) {
1151            _mapper->UseLookupTableScalarRangeOn();
1152            _mapper->SetLookupTable(_lut);
1153        }
1154        if (_contourMapper != NULL) {
1155            _contourMapper->UseLookupTableScalarRangeOn();
1156            _contourMapper->SetLookupTable(_lut);
1157        }
1158        _lut->DeepCopy(cmap->getLookupTable());
1159        switch (_colorMode) {
1160        case COLOR_CONSTANT:
1161        case COLOR_BY_SCALAR:
1162            _lut->SetRange(_dataRange);
1163            break;
1164        case COLOR_BY_VECTOR_MAGNITUDE:
1165            _lut->SetRange(_vectorMagnitudeRange);
1166            break;
1167        case COLOR_BY_VECTOR_X:
1168            _lut->SetRange(_vectorComponentRange[0]);
1169            break;
1170        case COLOR_BY_VECTOR_Y:
1171            _lut->SetRange(_vectorComponentRange[1]);
1172            break;
1173        case COLOR_BY_VECTOR_Z:
1174            _lut->SetRange(_vectorComponentRange[2]);
1175            break;
1176        default:
1177            break;
1178        }
1179    } else {
1180        double range[2];
1181        _lut->GetTableRange(range);
1182        _lut->DeepCopy(cmap->getLookupTable());
1183        _lut->SetRange(range);
1184        _lut->Modified();
1185    }
1186
1187    switch (_colorMode) {
1188    case COLOR_BY_VECTOR_MAGNITUDE:
1189        _lut->SetVectorModeToMagnitude();
1190        break;
1191    case COLOR_BY_VECTOR_X:
1192        _lut->SetVectorModeToComponent();
1193        _lut->SetVectorComponent(0);
1194        break;
1195    case COLOR_BY_VECTOR_Y:
1196        _lut->SetVectorModeToComponent();
1197        _lut->SetVectorComponent(1);
1198        break;
1199    case COLOR_BY_VECTOR_Z:
1200        _lut->SetVectorModeToComponent();
1201        _lut->SetVectorComponent(2);
1202        break;
1203    default:
1204         break;
1205    }
1206}
1207
1208/**
1209 * \brief Specify number of evenly spaced contour lines to render
1210 *
1211 * Will override any existing contours
1212 */
1213void HeightMap::setNumContours(int numContours)
1214{
1215    _contours.clear();
1216    _numContours = numContours;
1217
1218    update();
1219}
1220
1221/**
1222 * \brief Specify an explicit list of contour isovalues
1223 *
1224 * Will override any existing contours
1225 */
1226void HeightMap::setContourList(const std::vector<double>& contours)
1227{
1228    _contours = contours;
1229    _numContours = (int)_contours.size();
1230
1231    update();
1232}
1233
1234/**
1235 * \brief Get the number of contours
1236 */
1237int HeightMap::getNumContours() const
1238{
1239    return _numContours;
1240}
1241
1242/**
1243 * \brief Get the contour list (may be empty if number of contours
1244 * was specified in place of a list)
1245 */
1246const std::vector<double>& HeightMap::getContourList() const
1247{
1248    return _contours;
1249}
1250
1251/**
1252 * \brief Turn on/off lighting of this object
1253 */
1254void HeightMap::setLighting(bool state)
1255{
1256    _lighting = state;
1257    if (_dsActor != NULL)
1258        _dsActor->GetProperty()->SetLighting((state ? 1 : 0));
1259}
1260
1261/**
1262 * \brief Turn on/off rendering of mesh edges
1263 */
1264void HeightMap::setEdgeVisibility(bool state)
1265{
1266    if (_dsActor != NULL) {
1267        _dsActor->GetProperty()->SetEdgeVisibility((state ? 1 : 0));
1268    }
1269}
1270
1271/**
1272 * \brief Turn on/off rendering of colormaped surface
1273 */
1274void HeightMap::setContourSurfaceVisibility(bool state)
1275{
1276    if (_dsActor != NULL) {
1277        _dsActor->SetVisibility((state ? 1 : 0));
1278    }
1279}
1280
1281/**
1282 * \brief Turn on/off rendering of contour isolines
1283 */
1284void HeightMap::setContourLineVisibility(bool state)
1285{
1286    if (_contourActor != NULL) {
1287        _contourActor->SetVisibility((state ? 1 : 0));
1288    }
1289}
1290
1291/**
1292 * \brief Set RGB color of mesh
1293 */
1294void HeightMap::setColor(float color[3])
1295{
1296    _color[0] = color[0];
1297    _color[1] = color[1];
1298    _color[2] = color[2];
1299    if (_dsActor != NULL)
1300        _dsActor->GetProperty()->SetColor(_color[0], _color[1], _color[2]);
1301}
1302
1303/**
1304 * \brief Set RGB color of polygon edges
1305 */
1306void HeightMap::setEdgeColor(float color[3])
1307{
1308    _edgeColor[0] = color[0];
1309    _edgeColor[1] = color[1];
1310    _edgeColor[2] = color[2];
1311    if (_dsActor != NULL)
1312        _dsActor->GetProperty()->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
1313}
1314
1315/**
1316 * \brief Set pixel width of polygon edges (may be a no-op)
1317 */
1318void HeightMap::setEdgeWidth(float edgeWidth)
1319{
1320    _edgeWidth = edgeWidth;
1321    if (_dsActor != NULL)
1322        _dsActor->GetProperty()->SetLineWidth(_edgeWidth);
1323}
1324
1325/**
1326 * \brief Set RGB color of contour isolines
1327 */
1328void HeightMap::setContourEdgeColor(float color[3])
1329{
1330    _contourEdgeColor[0] = color[0];
1331    _contourEdgeColor[1] = color[1];
1332    _contourEdgeColor[2] = color[2];
1333    if (_contourActor != NULL) {
1334        _contourActor->GetProperty()->SetColor(_contourEdgeColor[0],
1335                                               _contourEdgeColor[1],
1336                                               _contourEdgeColor[2]);
1337        _contourActor->GetProperty()->SetEdgeColor(_contourEdgeColor[0],
1338                                                   _contourEdgeColor[1],
1339                                                   _contourEdgeColor[2]);
1340    }
1341}
1342
1343/**
1344 * \brief Set pixel width of contour isolines (may be a no-op)
1345 */
1346void HeightMap::setContourEdgeWidth(float edgeWidth)
1347{
1348    _contourEdgeWidth = edgeWidth;
1349    if (_contourActor != NULL)
1350        _contourActor->GetProperty()->SetLineWidth(_contourEdgeWidth);
1351}
1352
1353/**
1354 * \brief Set a group of world coordinate planes to clip rendering
1355 *
1356 * Passing NULL for planes will remove all cliping planes
1357 */
1358void HeightMap::setClippingPlanes(vtkPlaneCollection *planes)
1359{
1360    if (_mapper != NULL) {
1361        _mapper->SetClippingPlanes(planes);
1362    }
1363    if (_contourMapper != NULL) {
1364        _contourMapper->SetClippingPlanes(planes);
1365    }
1366}
1367
Note: See TracBrowser for help on using the repository browser.